source: soft/giet_vm/applications/fft/fft.c @ 812

Last change on this file since 812 was 812, checked in by alain, 8 years ago

Introduce application FFT.

File size: 35.3 KB
Line 
1/*************************************************************************/
2/*                                                                       */
3/*  Copyright (c) 1994 Stanford University                               */
4/*                                                                       */
5/*  All rights reserved.                                                 */
6/*                                                                       */
7/*  Permission is given to use, copy, and modify this software for any   */
8/*  non-commercial purpose as long as this copyright notice is not       */
9/*  removed.  All other uses, including redistribution in whole or in    */
10/*  part, are forbidden without prior written permission.                */
11/*                                                                       */
12/*  This software is provided with absolutely no warranty and no         */
13/*  support.                                                             */
14/*                                                                       */
15/*************************************************************************/
16
17///////////////////////////////////////////////////////////////////////////
18// This is the port of the SPLASH FFT benchmark on the GIET_VM OS.
19// It has been done by Hao LIU (04/12/2014), and Alain Greiner (14/03/2016).
20//
21// This application performs the 1D fast Fourier transfom for an array
22// of N complex points, using the Cooley-Tuckey FFT method.
23// The N data points are seen as a 2D array (rootN rows * rootN columns).
24// Each thread handle (rootN / nthreads) rows. The N input data points
25// be initialised in three different modes:
26// - CONSTANT : all data points have the same [1,0] value
27// - COSIN    : data point n has [cos(n/N) , sin(n/N)] values
28// - RANDOM   : data points have pseudo random values
29//
30// This application uses 4 shared data arrays, that are distributed
31// in all clusters (one sub-buffer per cluster):
32// - data[N] contains N input data points, with 2 double per point.
33// - trans[N] contains N intermediate data points, 2 double per point.
34// - umain[rootN] contains rootN coefs required for a rootN points FFT.
35// - twid[N] contains N coefs : exp(2*pi*i*j/N) / i and j in [0,rootN-1].
36// For data, trans, twid, each sub-buffer contains (N/nclusters) points.
37// For umain, each sub-buffer contains (rootN/nclusters) points.
38//
39// The two main parameters are the following:     
40//  - M : N = 2**M = number of data points / M must be an even number.
41//  - T : nthreads = nprocs defined by the hardware / must be power of 2.
42//
43// Several configuration parameters can be defined below:
44//  - VERBOSE : Print out complex data points arrays.
45//  - CHECK : Perform both FFT and inverse FFT to check output/input.
46//  - DEBUG : Display intermediate results
47///////////////////////////////////////////////////////////////////////////
48
49#include "stdio.h"
50#include "stdlib.h"
51#include "user_barrier.h"
52#include "malloc.h"
53#include "math.h"
54
55// constants
56
57#define PI                      3.14159265359
58#define PAGE_SIZE               4096
59#define CLUSTERS_MAX            256
60#define THREADS_MAX             1024
61#define RANDOM                  0
62#define COSIN                   1
63#define CONSTANT                2
64
65// parameters
66
67#define DEFAULT_M               16
68#define VERBOSE                 0
69#define CHECK                   0
70#define DEBUG                   0
71#define DETAILED_DEBUG          0
72#define MODE                    RANDOM
73
74// macro to swap two variables
75#define SWAP(a,b) { double tmp; tmp = a; a = b; b = tmp; }
76
77// global variables
78unsigned int   x_size;                     // number of clusters per row in the mesh
79unsigned int   y_size;                     // number of clusters per column in the mesh
80unsigned int   nprocs;                     // number of processors per cluster
81long           nthreads;                   // total number of threads (one thread per proc)
82long           nclusters;                  // total number of clusters
83long           M = DEFAULT_M;              // log2(number of points)
84long           N;                          // number of points (N = 2^M)         
85long           rootN;                      // rootN = 2^M/2   
86long           rows_per_proc;              // number of data "rows" handled by a single thread
87long           points_per_cluster;         // number of data points per cluster
88
89// arrays of pointers on distributed buffers (one sub-buffer per cluster)
90double *       data[CLUSTERS_MAX];         // original time-domain data
91double *       trans[CLUSTERS_MAX];        // used as auxiliary space for transpose
92double *       bloup[CLUSTERS_MAX];        // used as auxiliary space for DFT
93double *       umain[CLUSTERS_MAX];        // roots of unity used fo rootN points FFT   
94double *       twid[CLUSTERS_MAX];         // twiddle factor : exp(-2iPI*k*n/N)
95
96// instrumentation counters
97long           parallel_time[THREADS_MAX]; // total computation time (per thread)
98long           sync_time[THREADS_MAX];     // cumulative waiting time in barriers (per thread)
99long           init_time;                  // initialisation time (in main)
100
101// lock protecting shared TTY
102user_lock_t  tty_lock;
103
104// synchronisation barrier (all threads)
105giet_sqt_barrier_t barrier;
106
107// functions declaration
108void slave();
109
110double CheckSum(double ** x);
111
112void InitX(double ** x , unsigned int mode);
113
114void InitU(double ** u);
115
116void InitT(double ** u);
117
118long BitReverse( long k );
119
120void FFT1D( long direction , double ** x , double ** tmp , double * upriv, 
121            double ** twid , long MyNum , long MyFirst , long MyLast );
122
123void TwiddleOneCol(long direction, long j, double ** u, double ** x, long offset_x );
124
125void Scale( double **x, long offset_x );
126
127void Transpose( double ** src, double ** dest, long MyFirst, long MyLast );
128
129void Copy( double ** src, double ** dest, long MyFirst , long MyLast );
130
131void Reverse( double ** x, long offset_x );
132
133void FFT1DOnce( long direction , double * u , double ** x , long offset_x );
134
135void PrintArray( double ** x , long size );
136
137void SimpleDft( long direction , long size , double ** src , long src_offset ,
138                double ** dst , long dst_offset );
139
140///////////////////////////////////////////////////////////////////
141// This function execute the sequencial initialisation
142// launch the parallel execution, and makes the instrumentation.
143///////////////////////////////////////////////////////////////////
144__attribute__ ((constructor)) void main()
145{
146    long    x;             // index for cluster X coordinate
147    long    y;             // index for cluster Y coordinate
148    long    i;             // cluster index = y * x_size + x
149    long    n;             // thread index
150    long    start_time;    // starting cycle for main
151    double  ck1;           // for input/output checking
152    double  ck3;           // for input/output checking
153
154    // FFT application starts
155    start_time = giet_proctime();
156
157    // allocate shared TTY
158    giet_tty_alloc( 1 );     
159    lock_init( &tty_lock);
160
161    // get platform parameters to compute nthreads & nclusters
162    giet_procs_number( &x_size , &y_size , &nprocs );
163
164    giet_pthread_assert( ((nprocs == 1) || (nprocs == 2) || (nprocs == 4)),
165                         "[FFT ERROR] number of procs per cluster must be 1, 2 or 4");
166
167    giet_pthread_assert( ((x_size == 1) || (x_size == 2) || (x_size == 4) || 
168                          (x_size == 8) || (x_size == 16)),
169                         "[FFT ERROR] x_size must be 1,2,4,8,16");
170
171    giet_pthread_assert( ((y_size == 1) || (y_size == 2) || (y_size == 4) || 
172                          (y_size == 8) || (y_size == 16)),
173                         "[FFT ERROR] y_size must be 1,2,4,8,16");
174
175    nthreads  = x_size * y_size * nprocs;
176    nclusters = x_size * y_size;
177
178    // compute various constants depending on N and T
179    N                  = 1 << M;
180    rootN              = 1 << (M / 2);
181    rows_per_proc      = rootN / nthreads;
182    points_per_cluster = N / nclusters;
183 
184    giet_pthread_assert( (rootN >= nthreads) ,
185                         "[FFT ERROR] sqrt(N) must be larger than nthreads");
186
187    giet_tty_printf("\n*** Entering FFT  ***\n"
188                    " - %d Complex points\n"
189                    " - %d Thread(s)\n",
190                    N , nthreads );
191
192    // distributed heap initialisation
193    for ( x = 0 ; x < x_size ; x++ ) 
194    {
195        for ( y = 0 ; y < y_size ; y++ ) 
196        {
197            heap_init( x , y );
198        }
199    }
200
201    // allocate memory for the distributed data[i], trans[i], umain[i], twid[i] buffers
202    // the index (i) is a continuous cluster index
203    long data_size   = (N / nclusters) * 2 * sizeof(double);
204    long coefs_size  = (rootN / nclusters) * 2 * sizeof(double);
205    for (x = 0 ; x < x_size ; x++)
206    {
207        for (y = 0 ; y < y_size ; y++)
208        {
209            i         = x * y_size + y;
210            data[i]   = (double *)remote_malloc( data_size  , x , y ); 
211            trans[i]  = (double *)remote_malloc( data_size  , x , y ); 
212            bloup[i]  = (double *)remote_malloc( data_size  , x , y ); 
213            umain[i]  = (double *)remote_malloc( coefs_size , x , y ); 
214            twid[i]   = (double *)remote_malloc( data_size  , x , y ); 
215        }
216    }
217
218    // arrays initialisation
219    InitX( data , MODE ); 
220    InitU( umain ); 
221    InitT( twid );
222
223    if ( CHECK ) 
224    {
225        ck1 = CheckSum( data );
226    }
227
228    if ( VERBOSE ) 
229    {
230        giet_tty_printf("Data values / base = %x\n", &data[0][0] );
231        PrintArray( data , N );
232
233        giet_tty_printf("Twiddle values / base = %x\n", &twid[0][0] );
234        PrintArray( twid , N );
235
236        SimpleDft( 1 , N , data , 0 , bloup , 0 );
237        giet_tty_printf("Expected results / base = %x\n", &bloup[0][0] );
238        PrintArray( bloup , N );
239    }
240
241    // allocate thread[] array : thread identidiers defined by the kernel
242    pthread_t* thread = malloc( nthreads * sizeof(pthread_t) );
243
244    // allocate tid[] array : continuous thred index defined by the user
245    long* tid = malloc( nthreads * sizeof(long) );
246
247    // initialise distributed barrier
248    sqt_barrier_init( &barrier, x_size , y_size , nprocs );
249
250    // Initialisation completed
251    giet_tty_printf("\n*** Starting Parallel FFT  ***\n");
252
253    // register initialisation time
254    init_time = giet_proctime() - start_time;
255
256    // launch other threads to execute the slave() function
257    for ( n = 1 ; n < nthreads ; n++ )
258    {
259        tid[n] = n;
260        if ( giet_pthread_create( &thread[n],
261                                  NULL,                  // no attribute
262                                  &slave,
263                                  &tid[n] ) )            // pointer on thread index
264        {
265            giet_tty_printf("\n[FFT ERROR] creating thread %d\n", tid[n] );
266            giet_pthread_exit( NULL );
267        }
268    }
269
270    // main execute itself the slave() function
271    tid[0] = 0;
272    slave( &tid[0] );
273
274    // wait other threads completion
275    for ( n = 1 ; n < nthreads ; n++ )
276    {
277        if ( giet_pthread_join( thread[n], NULL ) )
278        {
279            giet_tty_printf("\n[FFT ERROR] joining thread %d\n", tid[n] );
280            giet_pthread_exit( NULL );
281        }
282    }
283
284    if ( VERBOSE ) 
285    {
286        giet_tty_printf("Data values after FFT:\n");
287        PrintArray( data , N );
288    }
289
290    if ( CHECK )
291    { 
292        ck3 = CheckSum( data );
293        giet_tty_printf("\n*** Results ***\n");
294        giet_tty_printf("Checksum difference is %f (%f, %f)\n", ck1 - ck3, ck1, ck3);
295        if (fabs(ck1 - ck3) < 0.001)  giet_tty_printf("Results OK.\n");
296        else                          giet_tty_printf("Results KO.\n");
297    }
298
299    // instrumentation (display and save on disk)
300    char string[256];
301
302    snprintf( string , 256 , "/home/fft_%d_%d_%d_%d",
303              x_size , y_size , nprocs , N );
304
305    // open instrumentation file
306    long fd = giet_fat_open( string , O_CREAT );
307    if ( fd < 0 ) 
308    { 
309        giet_tty_printf("\n[FFT ERROR] cannot open instrumentation file %s\n", string );
310        giet_pthread_exit( NULL );
311    }
312
313    snprintf( string , 256 , "\n*** FFT instrumentation : (%dx%dx%d) procs / %d data points ---\n",
314                             x_size, y_size, nprocs , N ); 
315
316    giet_tty_printf( "%s" , string );
317    giet_fat_fprintf( fd , "%s" , string );
318
319    long min_para = parallel_time[0];
320    long max_para = parallel_time[0];
321    long min_sync = sync_time[0];
322    long max_sync = sync_time[0];
323
324    for (i = 1 ; i < nthreads ; i++) 
325    {
326        if (parallel_time[i] > max_para)  max_para = parallel_time[i];
327        if (parallel_time[i] < min_para)  min_para = parallel_time[i];
328        if (sync_time[i]     > max_sync)  max_sync = sync_time[i];
329        if (sync_time[i]     < min_sync)  min_sync = sync_time[i];
330    }
331
332    snprintf( string , 256 , "\n      Init       Parallel   Barrier\n"
333                             "MIN : %d  |  %d  |  %d   (cycles)\n" 
334                             "MAX : %d  |  %d  |  %d   (cycles)\n",
335                             (int)init_time, (int)min_para, (int)min_sync,
336                             (int)init_time, (int)max_para, (int)max_sync );
337
338    giet_tty_printf("%s" , string );
339    giet_fat_fprintf( fd , "%s" , string );
340
341    // close instrumentation file and exit
342    giet_fat_close( fd );
343
344    giet_pthread_exit("main() completed");
345
346    return ;
347} // end main()
348
349///////////////////////////////////////////////////////////////
350// This function is executed in parallel by all threads.
351///////////////////////////////////////////////////////////////
352__attribute__ ((constructor)) void slave( long* tid ) 
353{
354    long     i;
355    long     MyNum = *tid;
356    long     MyFirst;         // index first row allocated to thread
357    long     MyLast;          // index last row allocated to thread
358    double * upriv;
359    long     c_id;
360    long     c_offset;
361    long     parallel_start;
362    long     barrier_start;
363
364    // initialise instrumentation
365    parallel_start   = giet_proctime();
366    sync_time[MyNum] = 0;
367
368    // allocate and initialise local array upriv[]
369    // that is a local copy of the rootN coefs defined in umain[]
370    upriv = (double *)malloc(2 * (rootN - 1) * sizeof(double)); 
371    for ( i = 0 ; i < (rootN - 1) ; i++) 
372    {
373        c_id     = i / (rootN / nclusters);
374        c_offset = i % (rootN / nclusters);
375        upriv[2*i]   = umain[c_id][2*c_offset];
376        upriv[2*i+1] = umain[c_id][2*c_offset+1];
377    }
378
379    // compute first and last rows handled by the thread
380    MyFirst = rootN * MyNum / nthreads;
381    MyLast  = rootN * (MyNum + 1) / nthreads;
382
383    // perform forward FFT
384    FFT1D( 1 , data , trans , upriv , twid , MyNum , MyFirst , MyLast );
385
386    // BARRIER
387    barrier_start = giet_proctime();
388    sqt_barrier_wait( &barrier );
389    sync_time[MyNum] = giet_proctime() - barrier_start;
390
391    // perform backward FFT if required
392    if ( CHECK ) 
393    {
394        // BARRIER
395        barrier_start = giet_proctime();
396        sqt_barrier_wait( &barrier );
397        sync_time[MyNum] = giet_proctime() - barrier_start;
398
399        FFT1D( -1 , data , trans , upriv , twid , MyNum , MyFirst , MyLast );
400    }
401
402    // register computation time
403    parallel_time[MyNum] = giet_proctime() - parallel_start;
404
405    // exit if MyNum != 0
406    if( MyNum ) giet_pthread_exit("slave() completed\n");
407
408}  // end slave()
409
410////////////////////////////////////////////////////////////////////////////////////////
411// This function makes the DFT from the src[nclusters][points_per_cluster] distributed
412// buffer, to the dst[nclusters][points_per_cluster] distributed buffer.
413////////////////////////////////////////////////////////////////////////////////////////
414void SimpleDft( long      direction,
415                long      size,           // number of points
416                double ** src,            // source distributed buffer
417                long      src_offset,     // offset in source array
418                double ** dst,            // destination distributed buffer
419                long      dst_offset )    // offset in destination array
420{
421    long    n , k;
422    double  phi;            // 2*PI*n*k/N
423    double  u_r;            // cos( phi )
424    double  u_c;            // sin( phi )
425    double  d_r;            // Re(data[n])
426    double  d_c;            // Im(data[n])
427    double  accu_r;         // Re(accu)
428    double  accu_c;         // Im(accu)
429    long    c_id;           // distributed buffer cluster index
430    long    c_offset;       // offset in distributed buffer
431
432    for ( k = 0 ; k < size ; k++ )       // loop on the output data points
433    {
434        // initialise accu
435        accu_r = 0;
436        accu_c = 0;
437
438        for ( n = 0 ; n < size ; n++ )   // loop on the input data points
439        {
440            // compute coef
441            phi = (double)(2*PI*n*k) / size;
442            u_r =  cos( phi );
443            u_c = -sin( phi ) * direction;
444
445            // get input data point
446            c_id     = (src_offset + n) / (points_per_cluster);
447            c_offset = (src_offset + n) % (points_per_cluster);
448            d_r      = data[c_id][2*c_offset];
449            d_c      = data[c_id][2*c_offset+1];
450
451            // increment accu
452            accu_r += ((u_r*d_r) - (u_c*d_c));
453            accu_c += ((u_r*d_c) + (u_c*d_r));
454        }
455
456        // scale for inverse DFT
457        if ( direction == -1 )
458        {
459            accu_r /= size;
460            accu_c /= size;
461        }
462
463        // set output data point
464        c_id     = (dst_offset + k) / (points_per_cluster);
465        c_offset = (dst_offset + k) % (points_per_cluster);
466        dst[c_id][2*c_offset]   = accu_r;
467        dst[c_id][2*c_offset+1] = accu_c;
468    }
469
470}  // end SimpleDft()
471
472////////////////////////////
473double CheckSum(double ** x) 
474{
475    long i , j;
476    double       cks;
477    long c_id;
478    long c_offset;
479
480    cks = 0.0;
481    for (j = 0; j < rootN ; j++) 
482    {
483        for (i = 0; i < rootN ; i++) 
484        {
485            c_id      = (rootN * j + i) / (points_per_cluster);
486            c_offset  = (rootN * j + i) % (points_per_cluster);
487
488            cks += data[c_id][2*c_offset] + data[c_id][2*c_offset+1];
489        }
490    }
491    return(cks);
492}
493
494
495///////////////////////
496void InitX(double ** x,
497           unsigned int mode ) 
498{
499    long    i , j;
500    long    c_id;
501    long    c_offset;
502    long    index;
503
504    for ( j = 0 ; j < rootN ; j++ )      // loop on row index
505    { 
506        for ( i = 0 ; i < rootN ; i++ )  // loop on point in a row
507        { 
508            index     = j * rootN + i;
509            c_id      = index / (points_per_cluster);
510            c_offset  = index % (points_per_cluster);
511
512            // complex input signal is random
513            if ( mode == RANDOM )               
514            {
515                data[c_id][2*c_offset]   = ( (double)giet_rand() ) / 65536;
516                data[c_id][2*c_offset+1] = ( (double)giet_rand() ) / 65536;
517            }
518           
519
520            // complex input signal is cos(n/N) / sin(n/N)
521            if ( mode == COSIN )               
522            {
523                double phi = (double)( 2 * PI * index) / N;
524                data[c_id][2*c_offset]   = cos( phi );
525                data[c_id][2*c_offset+1] = sin( phi );
526            }
527
528            // complex input signal is constant
529            if ( mode == CONSTANT )               
530            {
531                data[c_id][2*c_offset]   = 1.0;
532                data[c_id][2*c_offset+1] = 0.0;
533            }
534        }
535    }
536}
537
538/////////////////////////
539void InitU( double ** u ) 
540{
541    long    q; 
542    long    j; 
543    long    base; 
544    long    n1;
545    long    c_id;
546    long    c_offset;
547    double  phi;
548    long    stop = 0;
549
550    for (q = 0 ; ((1 << q) < N) && (stop == 0) ; q++) 
551    { 
552        n1 = 1 << q;
553        base = n1 - 1;
554        for (j = 0; (j < n1) && (stop == 0) ; j++) 
555        {
556            if (base + j > rootN - 1) return;
557
558            c_id      = (base + j) / (rootN / nclusters);
559            c_offset  = (base + j) % (rootN / nclusters);
560            phi = (double)(2.0 * PI * j) / (2 * n1);
561            u[c_id][2*c_offset]   = cos( phi );
562            u[c_id][2*c_offset+1] = -sin( phi );
563        }
564    }
565}
566
567//////////////////////////
568void InitT( double ** u )
569{
570    long    i, j;
571    long    index;
572    long    c_id;
573    long    c_offset;
574    double  phi;
575
576    for ( j = 0 ; j < rootN ; j++ )      // loop on row index
577    { 
578        for ( i = 0 ; i < rootN ; i++ )  // loop on points in a row
579        { 
580            index     = j * rootN + i;
581            c_id      = index / (points_per_cluster);
582            c_offset  = index % (points_per_cluster);
583
584            phi = (double)(2.0 * PI * i * j) / N;
585            u[c_id][2*c_offset]   = cos( phi );
586            u[c_id][2*c_offset+1] = -sin( phi );
587        }
588    }
589}
590
591////////////////////////////////////////////////////////////////////////////////////////
592// This function returns an index value that is the bit reverse of the input value.
593////////////////////////////////////////////////////////////////////////////////////////
594long BitReverse( long k ) 
595{
596    long i; 
597    long j; 
598    long tmp;
599
600    j = 0;
601    tmp = k;
602    for (i = 0; i < M/2 ; i++) 
603    {
604        j = 2 * j + (tmp & 0x1);
605        tmp = tmp >> 1;
606    }
607    return j;
608}
609
610////////////////////////////////////////////////////////////////////////////////////////
611// This function perform the in place (direct or inverse) FFT on the N data points
612// contained in the distributed buffers x[nclusters][points_per_cluster].
613// It handles the (N) points 1D array as a (rootN*rootN) points 2D array.
614// 1) it transpose (rootN/nthreads ) rows from x to tmp.
615// 2) it make (rootN/nthreads) FFT on the tmp rows and apply the twiddle factor.
616// 3) it transpose (rootN/nthreads) columns from tmp to x.
617// 4) it make (rootN/nthreads) FFT on the x rows.
618// It calls the FFT1DOnce() 2*(rootN/nthreads) times to perform the in place FFT
619// on the rootN points contained in a row.
620////////////////////////////////////////////////////////////////////////////////////////
621void FFT1D( long       direction,       // direct : 1 / inverse : -1
622            double **  x,               // input & output distributed data points array
623            double **  tmp,             // auxiliary distributed data points array
624            double *   upriv,           // local array containing coefs for rootN FFT
625            double **  twid,            // distributed arrays containing N twiddle factors
626            long       MyNum,
627            long       MyFirst, 
628            long       MyLast )
629{
630    long j;
631    long barrier_start;
632
633    // transpose (rootN/nthreads) rows from x to tmp
634    Transpose( x , tmp , MyFirst , MyLast );
635
636#if DEBUG
637giet_tty_printf("\n@@@ tmp after first transpose\n");
638PrintArray( tmp , N );
639#endif
640
641    // BARRIER
642    barrier_start = giet_proctime();
643    sqt_barrier_wait( &barrier );
644    sync_time[MyNum] = giet_proctime() - barrier_start;
645
646    // do FFTs on rows of tmp (i.e. columns of x) and apply twiddle factor
647    for (j = MyFirst; j < MyLast; j++) 
648    {
649        FFT1DOnce( direction , upriv , tmp , j * rootN );
650        TwiddleOneCol( direction , j , twid , tmp , j * rootN );
651    } 
652
653#if DEBUG
654giet_tty_printf("\n@@@ tmp after columns FFT + twiddle \n");
655PrintArray( tmp , N );
656#endif
657
658    // BARRIER
659    barrier_start = giet_proctime();
660    sqt_barrier_wait( &barrier );
661    sync_time[MyNum] = giet_proctime() - barrier_start;
662
663    // transpose tmp to x
664    Transpose( tmp , x , MyFirst , MyLast );
665
666#if DEBUG
667giet_tty_printf("\n@@@ x after second transpose \n");
668PrintArray( x , N );
669#endif
670
671    // BARRIER
672    barrier_start = giet_proctime();
673    sqt_barrier_wait( &barrier );
674    sync_time[MyNum] = giet_proctime() - barrier_start;
675
676    // do FFTs on rows of x and apply the scaling factor
677    for (j = MyFirst; j < MyLast; j++) 
678    {
679        FFT1DOnce( direction , upriv , x , j * rootN );
680        if (direction == -1) Scale( x , j * rootN );
681    }
682
683#if DEBUG
684giet_tty_printf("\n@@@ x after rows FFT + scaling \n");
685PrintArray( x , N );
686#endif
687
688    // BARRIER
689    barrier_start = giet_proctime();
690    sqt_barrier_wait( &barrier );
691    sync_time[MyNum] = giet_proctime() - barrier_start;
692
693    // transpose x to tmp
694    Transpose( x , tmp , MyFirst , MyLast );
695
696#if DEBUG
697giet_tty_printf("\n@@@ tmp after third transpose \n");
698PrintArray( tmp , N );
699#endif
700
701    // BARRIER
702    barrier_start = giet_proctime();
703    sqt_barrier_wait( &barrier );
704    sync_time[MyNum] = giet_proctime() - barrier_start;
705
706    // copy tmp to x
707    Copy( tmp , x , MyFirst , MyLast );
708
709#if DEBUG
710giet_tty_printf("\n@@@ x after final copy \n");
711PrintArray( x , N );
712#endif
713
714
715}  // end FFT1D()
716
717/////////////////////////////////////////////////////////////////////////////////////
718// This function multiply all points contained in a row (rootN points) of the
719// x[] array by the corresponding twiddle factor, contained in the u[] array.
720/////////////////////////////////////////////////////////////////////////////////////
721void TwiddleOneCol( long      direction, 
722                    long      j,              // y coordinate in 2D view of coef array
723                    double ** u,              // coef array base address
724                    double ** x,              // data array base address
725                    long      offset_x )      // first point in N points data array
726{
727    long i;
728    double       omega_r; 
729    double       omega_c; 
730    double       x_r; 
731    double       x_c;
732    long         c_id;
733    long         c_offset;
734
735    for (i = 0; i < rootN ; i++)  // loop on the rootN points
736    {
737        // get coef
738        c_id      = (j * rootN + i) / (points_per_cluster);
739        c_offset  = (j * rootN + i) % (points_per_cluster);
740        omega_r = u[c_id][2*c_offset];
741        omega_c = direction * u[c_id][2*c_offset+1];
742
743        // access data
744        c_id      = (offset_x + i) / (points_per_cluster);
745        c_offset  = (offset_x + i) % (points_per_cluster);   
746        x_r = x[c_id][2*c_offset]; 
747        x_c = x[c_id][2*c_offset+1];
748
749        x[c_id][2*c_offset]   = omega_r*x_r - omega_c * x_c;
750        x[c_id][2*c_offset+1] = omega_r*x_c + omega_c * x_r;
751    }
752}  // end TwiddleOneCol()
753
754////////////////////////
755void Scale( double ** x,           // data array base address
756            long      offset_x )   // first point of the row to be scaled
757{
758    long i;
759    long c_id;
760    long c_offset;
761
762    for (i = 0; i < rootN ; i++) 
763    {
764        c_id      = (offset_x + i) / (points_per_cluster);
765        c_offset  = (offset_x + i) % (points_per_cluster);
766        data[c_id][2*c_offset]     /= N;
767        data[c_id][2*c_offset + 1] /= N;
768    }
769}
770
771////////////////////////////
772void Transpose( double ** src,      // source buffer (array of pointers)
773                double ** dest,     // destination buffer (array of pointers)
774                long      MyFirst,  // first row allocated to the thread
775                long      MyLast )  // last row allocated to the thread
776{
777    long row;               // row index
778    long point;             // data point index in a row
779
780    long index_src;         // absolute index in the source N points array
781    long c_id_src;          // cluster for the source buffer
782    long c_offset_src;      // offset in the source buffer
783
784    long index_dst;         // absolute index in the dest N points array
785    long c_id_dst;          // cluster for the dest buffer
786    long c_offset_dst;      // offset in the dest buffer
787
788   
789    // scan all data points allocated to the thread
790    // (between MyFirst row and MyLast row) from the source buffer
791    // and write these points to the destination buffer
792    for ( row = MyFirst ; row < MyLast ; row++ )       // loop on the rows
793    {
794        for ( point = 0 ; point < rootN ; point++ )    // loop on points in row
795        {
796            index_src    = row * rootN + point;
797            c_id_src     = index_src / (points_per_cluster);
798            c_offset_src = index_src % (points_per_cluster);
799
800            index_dst    = point * rootN + row;
801            c_id_dst     = index_dst / (points_per_cluster);
802            c_offset_dst = index_dst % (points_per_cluster);
803
804            dest[c_id_dst][2*c_offset_dst]   = src[c_id_src][2*c_offset_src];
805            dest[c_id_dst][2*c_offset_dst+1] = src[c_id_src][2*c_offset_src+1];
806        }
807    }
808}  // end Transpose()
809
810/////////////////////////
811void Copy( double ** src,      // source buffer (array of pointers)
812           double ** dest,     // destination buffer (array of pointers)
813           long      MyFirst,  // first row allocated to the thread
814           long      MyLast )  // last row allocated to the thread
815{
816    long row;                  // row index
817    long point;                // data point index in a row
818
819    long index;                // absolute index in the N points array
820    long c_id;                 // cluster index
821    long c_offset;             // offset in local buffer
822
823    // scan all data points allocated to the thread
824    for ( row = MyFirst ; row < MyLast ; row++ )       // loop on the rows
825    {
826        for ( point = 0 ; point < rootN ; point++ )    // loop on points in row
827        {
828            index    = row * rootN + point;
829            c_id     = index / (points_per_cluster);
830            c_offset = index % (points_per_cluster);
831
832            dest[c_id][2*c_offset]   = src[c_id][2*c_offset];
833            dest[c_id][2*c_offset+1] = src[c_id][2*c_offset+1];
834        }
835    }
836}  // end Copy()
837
838//////////////////////////
839void Reverse( double ** x, 
840              long      offset_x )
841{
842    long j, k;
843    long c_id_j;
844    long c_offset_j;
845    long c_id_k;
846    long c_offset_k;
847
848    for (k = 0 ; k < rootN ; k++) 
849    {
850        j = BitReverse( k );
851        if (j > k) 
852        {
853            c_id_j      = (offset_x + j) / (points_per_cluster);
854            c_offset_j  = (offset_x + j) % (points_per_cluster);
855            c_id_k      = (offset_x + k) / (points_per_cluster);
856            c_offset_k  = (offset_x + k) % (points_per_cluster);
857
858            SWAP(x[c_id_j][2*c_offset_j]  , x[c_id_k][2*c_offset_k]);
859            SWAP(x[c_id_j][2*c_offset_j+1], x[c_id_k][2*c_offset_k+1]);
860        }
861    }
862}
863
864/////////////////////////////////////////////////////////////////////////////
865// This function makes the in-place FFT on all points contained in a row
866// (i.e. rootN points) of the x[nclusters][points_per_cluster] array.
867/////////////////////////////////////////////////////////////////////////////
868void FFT1DOnce( long      direction,  // direct / inverse
869                double *  u,          // private coefs array
870                double ** x,          // array of pointers on distributed buffers
871                long      offset_x )  // absolute offset in the x array
872{
873    long     j;
874    long     k;
875    long     q;
876    long     L;
877    long     r;
878    long     Lstar;
879    double * u1; 
880
881    long     offset_x1;     // index first butterfly input
882    long     offset_x2;     // index second butterfly output
883
884    double   omega_r;       // real part butterfy coef
885    double   omega_c;       // complex part butterfly coef
886
887    double   tau_r;
888    double   tau_c;
889
890    double   d1_r;          // real part first butterfly input
891    double   d1_c;          // imag part first butterfly input
892    double   d2_r;          // real part second butterfly input
893    double   d2_c;          // imag part second butterfly input
894
895    long     c_id_1;        // cluster index for first butterfly input
896    long     c_offset_1;    // offset for first butterfly input
897    long     c_id_2;        // cluster index for second butterfly input
898    long     c_offset_2;    // offset for second butterfly input
899
900#if DETAILED_DEBUG
901unsigned int p;
902giet_tty_printf("\n@@@ FFT ROW data in / %d points / offset = %d\n",
903                rootN , offset_x );
904for ( p = 0 ; p < rootN ; p++ )
905{
906    long index    = offset_x + p;
907    long c_id     = index / (points_per_cluster);
908    long c_offset = index % (points_per_cluster);
909    giet_tty_printf("%f , %f | ", x[c_id][2*c_offset] , x[c_id][2*c_offset+1] );
910}
911giet_tty_printf("\n");
912#endif
913
914    // This makes the rootN input points reordering
915    Reverse( x , offset_x ); 
916
917#if DETAILED_DEBUG
918giet_tty_printf("\n@@@ FFT ROW data after reverse\n");
919for ( p = 0 ; p < rootN ; p++ )
920{
921    long index    = offset_x + p;
922    long c_id     = index / (points_per_cluster);
923    long c_offset = index % (points_per_cluster);
924    giet_tty_printf("%f , %f | ", x[c_id][2*c_offset] , x[c_id][2*c_offset+1] );
925}
926giet_tty_printf("\n");
927#endif
928
929    // This implements the multi-stages, in place Butterfly network
930    for (q = 1; q <= M/2 ; q++)     // loop on stages
931    {
932        L = 1 << q;       // number of points per subset for current stage
933        r = rootN / L;    // number of subsets
934        Lstar = L / 2;
935        u1 = &u[2 * (Lstar - 1)];
936        for (k = 0; k < r; k++)     // loop on the subsets
937        {
938            offset_x1  = offset_x + (k * L);            // index first point
939            offset_x2  = offset_x + (k * L + Lstar);    // index second point
940
941#if DETAILED_DEBUG
942giet_tty_printf("\n ### q = %d / k = %d / x1 = %d / x2 = %d\n",
943                 q , k , offset_x1 , offset_x2 );
944#endif
945            // makes all in-place butterfly(s) for subset
946            for (j = 0; j < Lstar; j++) 
947            {
948                // get coef
949                omega_r = u1[2*j];
950                omega_c = direction * u1[2*j+1];
951
952                // get d[x1] address and value
953                c_id_1      = (offset_x1 + j) / (points_per_cluster);
954                c_offset_1  = (offset_x1 + j) % (points_per_cluster);
955                d1_r        = x[c_id_1][2*c_offset_1];
956                d1_c        = x[c_id_1][2*c_offset_1+1];
957
958                // get d[x2] address and value
959                c_id_2      = (offset_x2 + j) / (points_per_cluster);
960                c_offset_2  = (offset_x2 + j) % (points_per_cluster);
961                d2_r        = x[c_id_2][2*c_offset_2];
962                d2_c        = x[c_id_2][2*c_offset_2+1];
963
964#if DETAILED_DEBUG
965giet_tty_printf("\n ### d1_in = (%f , %f) / d2_in = (%f , %f) / coef = (%f , %f)\n", 
966                d1_r , d1_c , d2_r , d2_c , omega_r , omega_c);
967#endif
968                // tau = omega * d[x2]
969                tau_r = omega_r * d2_r - omega_c * d2_c;
970                tau_c = omega_r * d2_c + omega_c * d2_r;
971
972                // set new value for d[x1] = d[x1] + omega * d[x2]
973                x[c_id_1][2*c_offset_1]   = d1_r + tau_r;
974                x[c_id_1][2*c_offset_1+1] = d1_c + tau_c;
975
976                // set new value for d[x2] = d[x1] - omega * d[x2]
977                x[c_id_2][2*c_offset_2]   = d1_r - tau_r;
978                x[c_id_2][2*c_offset_2+1] = d1_c - tau_c;
979
980#if DETAILED_DEBUG
981giet_tty_printf("\n ### d1_out = (%f , %f) / d2_out = (%f , %f)\n", 
982                d1_r + tau_r , d1_c + tau_c , d2_r - tau_r , d2_c - tau_c );
983#endif
984            }
985        }
986    }
987
988#if DETAILED_DEBUG
989giet_tty_printf("\n@@@ FFT ROW data out\n");
990for ( p = 0 ; p < rootN ; p++ )
991{
992    long index    = offset_x + p;
993    long c_id     = index / (points_per_cluster);
994    long c_offset = index % (points_per_cluster);
995    giet_tty_printf("%f , %f | ", x[c_id][2*c_offset] , x[c_id][2*c_offset+1] );
996}
997giet_tty_printf("\n");
998#endif
999
1000}  // end FFT1DOnce()
1001
1002//////////////////////////////////
1003void PrintArray( double ** array,
1004                 long      size ) 
1005{
1006    long  i;
1007    long  c_id;
1008    long  c_offset;
1009
1010    // float display
1011    for (i = 0; i < size ; i++) 
1012    {
1013        c_id      = i / (points_per_cluster);
1014        c_offset  = i % (points_per_cluster);
1015
1016        giet_tty_printf(" %f  %f |", 
1017                        array[c_id][2*c_offset],
1018                        array[c_id][2*c_offset+1]);
1019
1020        if ( (i+1) % 4 == 0)  giet_tty_printf("\n");
1021    }
1022    giet_tty_printf("\n");
1023}
1024
1025
1026// Local Variables:
1027// tab-width: 4
1028// c-basic-offset: 4
1029// c-file-offsets:((innamespace . 0)(inline-open . 0))
1030// indent-tabs-mode: nil
1031// End:
1032
1033// vim: filetype=cpp:expandtab:shiftwidth=4:tabstop=4:softtabstop=4
1034
Note: See TracBrowser for help on using the repository browser.