space.c 39.8 KB
Newer Older
Pedro Gonnet's avatar
Pedro Gonnet committed
1
/*******************************************************************************
2
 * This file is part of SWIFT.
Pedro Gonnet's avatar
Pedro Gonnet committed
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
 * Coypright (c) 2012 Pedro Gonnet (pedro.gonnet@durham.ac.uk)
 * 
 * This program is free software: you can redistribute it and/or modify
 * it under the terms of the GNU Lesser General Public License as published
 * by the Free Software Foundation, either version 3 of the License, or
 * (at your option) any later version.
 * 
 * This program is distributed in the hope that it will be useful,
 * but WITHOUT ANY WARRANTY; without even the implied warranty of
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 * GNU General Public License for more details.
 * 
 * You should have received a copy of the GNU Lesser General Public License
 * along with this program.  If not, see <http://www.gnu.org/licenses/>.
 * 
 ******************************************************************************/

/* Config parameters. */
#include "../config.h"

/* Some standard headers. */
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <pthread.h>
#include <float.h>
#include <limits.h>
#include <math.h>
31
#include <omp.h>
Pedro Gonnet's avatar
Pedro Gonnet committed
32

33
34
35
36
37
/* MPI headers. */
#ifdef WITH_MPI
    #include <mpi.h>
#endif

Pedro Gonnet's avatar
Pedro Gonnet committed
38
/* Local headers. */
39
#include "const.h"
Pedro Gonnet's avatar
Pedro Gonnet committed
40
#include "cycle.h"
41
#include "atomic.h"
Pedro Gonnet's avatar
Pedro Gonnet committed
42
#include "lock.h"
43
#include "task.h"
44
#include "kernel.h"
45
#include "part.h"
Pedro Gonnet's avatar
Pedro Gonnet committed
46
#include "space.h"
47
#include "multipole.h"
48
49
50
#include "cell.h"
#include "scheduler.h"
#include "engine.h"
Pedro Gonnet's avatar
Pedro Gonnet committed
51
#include "runner.h"
52
#include "error.h"
Pedro Gonnet's avatar
Pedro Gonnet committed
53
54
55

/* Split size. */
int space_splitsize = space_splitsize_default;
56
int space_subsize = space_subsize_default;
57
int space_maxsize = space_maxsize_default;
Pedro Gonnet's avatar
Pedro Gonnet committed
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88

/* Map shift vector to sortlist. */
const int sortlistID[27] = {
    /* ( -1 , -1 , -1 ) */   0 ,
    /* ( -1 , -1 ,  0 ) */   1 , 
    /* ( -1 , -1 ,  1 ) */   2 ,
    /* ( -1 ,  0 , -1 ) */   3 ,
    /* ( -1 ,  0 ,  0 ) */   4 , 
    /* ( -1 ,  0 ,  1 ) */   5 ,
    /* ( -1 ,  1 , -1 ) */   6 ,
    /* ( -1 ,  1 ,  0 ) */   7 , 
    /* ( -1 ,  1 ,  1 ) */   8 ,
    /* (  0 , -1 , -1 ) */   9 ,
    /* (  0 , -1 ,  0 ) */   10 , 
    /* (  0 , -1 ,  1 ) */   11 ,
    /* (  0 ,  0 , -1 ) */   12 ,
    /* (  0 ,  0 ,  0 ) */   0 , 
    /* (  0 ,  0 ,  1 ) */   12 ,
    /* (  0 ,  1 , -1 ) */   11 ,
    /* (  0 ,  1 ,  0 ) */   10 , 
    /* (  0 ,  1 ,  1 ) */   9 ,
    /* (  1 , -1 , -1 ) */   8 ,
    /* (  1 , -1 ,  0 ) */   7 , 
    /* (  1 , -1 ,  1 ) */   6 ,
    /* (  1 ,  0 , -1 ) */   5 ,
    /* (  1 ,  0 ,  0 ) */   4 , 
    /* (  1 ,  0 ,  1 ) */   3 ,
    /* (  1 ,  1 , -1 ) */   2 ,
    /* (  1 ,  1 ,  0 ) */   1 , 
    /* (  1 ,  1 ,  1 ) */   0 
    };
89
90
    
    
91
92
93
94
95
96
97
98
99
100
101
102
103
104
/**
 * @brief Get the shift-id of the given pair of cells, swapping them
 *      if need be.
 *
 * @param s The space
 * @param ci Pointer to first #cell.
 * @param cj Pointer second #cell.
 * @param shift Vector from ci to cj.
 *
 * @return The shift ID and set shift, may or may not swap ci and cj.
 */
 
int space_getsid ( struct space *s , struct cell **ci , struct cell **cj , double *shift ) {

105
    int k, sid = 0, periodic = s->periodic;
106
107
108
109
110
111
    struct cell *temp;
    double dx[3];

    /* Get the relative distance between the pairs, wrapping. */
    for ( k = 0 ; k < 3 ; k++ ) {
        dx[k] = (*cj)->loc[k] - (*ci)->loc[k];
112
        if ( periodic && dx[k] < -s->dim[k]/2 )
113
            shift[k] = s->dim[k];
114
        else if ( periodic && dx[k] > s->dim[k]/2 )
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
            shift[k] = -s->dim[k];
        else
            shift[k] = 0.0;
        dx[k] += shift[k];
        }
        
    /* Get the sorting index. */
    for ( k = 0 ; k < 3 ; k++ )
        sid = 3*sid + ( (dx[k] < 0.0) ? 0 : ( (dx[k] > 0.0) ? 2 : 1 ) );

    /* Switch the cells around? */
    if ( runner_flip[sid] ) {
        temp = *ci; *ci = *cj; *cj = temp;
        for ( k = 0 ; k < 3 ; k++ )
            shift[k] = -shift[k];
        }
    sid = sortlistID[sid];
    
    /* Return the sort ID. */
    return sid;

    }


139
/**
140
 * @brief Recursively dismantle a cell tree.
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
 *
 */
 
void space_rebuild_recycle ( struct space *s , struct cell *c ) {
    
    int k;
    
    if ( c->split )
        for ( k = 0 ; k < 8 ; k++ )
            if ( c->progeny[k] != NULL ) {
                space_rebuild_recycle( s , c->progeny[k] );
                space_recycle( s , c->progeny[k] );
                c->progeny[k] = NULL;
                }
    
    }
157
158
    
    
159
/**
160
 * @brief Re-build the cell grid.
161
 *
162
163
 * @param s The #space.
 * @param cell_max Maximum cell edge length.
164
165
 */
 
166
void space_regrid ( struct space *s , double cell_max ) {
167

168
    float h_max = s->cell_min / kernel_gamma / space_stretch, dmin;
169
    int i, j, k, cdim[3], nr_parts = s->nr_parts;
170
    struct cell *restrict c;
171
    // ticks tic;
172
173
    
    /* Run through the parts and get the current h_max. */
174
    // tic = getticks();
175
176
177
178
179
180
181
182
183
184
185
186
    if ( s->cells != NULL ) {
        for ( k = 0 ; k < s->nr_cells ; k++ ) {
            if ( s->cells[k].h_max > h_max )
                h_max = s->cells[k].h_max;
            }
        }
    else {
        for ( k = 0 ; k < nr_parts ; k++ ) {
            if ( s->parts[k].h > h_max )
                h_max = s->parts[k].h;
            }
        s->h_max = h_max;
187
        }
188
189
190
191
192
193
194
195
196
197
198
199
        
    /* If we are running in parallel, make sure everybody agrees on
       how large the largest cell should be. */
    #ifdef WITH_MPI
    {
      float buff;
      if ( MPI_Allreduce( &h_max , &buff , 1 , MPI_FLOAT , MPI_MAX , MPI_COMM_WORLD ) != MPI_SUCCESS )
          error( "Failed to aggreggate the rebuild flag accross nodes." );
      h_max = buff;
    }
    #endif
    message( "h_max is %.3e (cell_max=%.3e)." , h_max , cell_max );
200
201
202
    
    /* Get the new putative cell dimensions. */
    for ( k = 0 ; k < 3 ; k++ )
203
        cdim[k] = floor( s->dim[k] / fmax( h_max*kernel_gamma*space_stretch , cell_max ) );
204
        
205
206
207
208
    /* Check if we have enough cells for periodicity. */
    if ( s->periodic && (cdim[0] < 3 || cdim[1] < 3 || cdim[2] < 3) )
        error( "Must have at least 3 cells in each spatial dimension when periodicity is switched on." );
        
209
210
211
212
213
214
    /* In MPI-Land, we're not allowed to change the top-level cell size. */
    #ifdef WITH_MPI
        if ( cdim[0] < s->cdim[0] || cdim[1] < s->cdim[1] || cdim[2] < s->cdim[2] )
            error( "Root-level change of cell size not allowed." );
    #endif
        
215
    /* Do we need to re-build the upper-level cells? */
216
    // tic = getticks();
Pedro Gonnet's avatar
Pedro Gonnet committed
217
    if ( s->cells == NULL ||
218
         cdim[0] < s->cdim[0] || cdim[1] < s->cdim[1] || cdim[2] < s->cdim[2] ) {
219
220
221
    
        /* Free the old cells, if they were allocated. */
        if ( s->cells != NULL ) {
222
            for ( k = 0 ; k < s->nr_cells ; k++ ) {
223
                space_rebuild_recycle( s , &s->cells[k] );
224
225
226
                if ( s->cells[k].sort != NULL )
                    free( s->cells[k].sort );
                }
227
228
229
230
            free( s->cells );
            s->maxdepth = 0;
            }
            
231
        /* Set the new cell dimensions only if smaller. */
232
233
234
235
236
        for ( k = 0 ; k < 3 ; k++ ) {
            s->cdim[k] = cdim[k];
            s->h[k] = s->dim[k] / cdim[k];
            s->ih[k] = 1.0 / s->h[k];
            }
237
        dmin = fminf( s->h[0] , fminf( s->h[1] , s->h[2] ) );
238

239
        /* Allocate the highest level of cells. */
Pedro Gonnet's avatar
Pedro Gonnet committed
240
        s->tot_cells = s->nr_cells = cdim[0] * cdim[1] * cdim[2];
241
242
243
244
245
246
247
248
249
250
251
252
253
254
        if ( posix_memalign( (void *)&s->cells , 64 , s->nr_cells * sizeof(struct cell) ) != 0 )
            error( "Failed to allocate cells." );
        bzero( s->cells , s->nr_cells * sizeof(struct cell) );
        for ( k = 0 ; k < s->nr_cells ; k++ )
            if ( lock_init( &s->cells[k].lock ) != 0 )
                error( "Failed to init spinlock." );

        /* Set the cell location and sizes. */
        for ( i = 0 ; i < cdim[0] ; i++ )
            for ( j = 0 ; j < cdim[1] ; j++ )
                for ( k = 0 ; k < cdim[2] ; k++ ) {
                    c = &s->cells[ cell_getid( cdim , i , j , k ) ];
                    c->loc[0] = i*s->h[0]; c->loc[1] = j*s->h[1]; c->loc[2] = k*s->h[2];
                    c->h[0] = s->h[0]; c->h[1] = s->h[1]; c->h[2] = s->h[2];
255
                    c->dmin = dmin;
256
                    c->depth = 0;
Pedro Gonnet's avatar
Pedro Gonnet committed
257
                    c->count = 0;
258
                    c->gcount = 0;
259
                    c->super = c;
260
                    lock_init( &c->lock );
261
                    }
262
263
           
        /* Be verbose about the change. */         
264
        message( "set cell dimensions to [ %i %i %i ]." , cdim[0] , cdim[1] , cdim[2] ); fflush(stdout);
265
266
                    
        } /* re-build upper-level cells? */
267
    // message( "rebuilding upper-level cells took %.3f ms." , (double)(getticks() - tic) / CPU_TPS * 1000 );
268
        
Pedro Gonnet's avatar
Pedro Gonnet committed
269
270
271
272
273
274
275
276
277
    /* Otherwise, just clean up the cells. */
    else {
    
        /* Free the old cells, if they were allocated. */
        for ( k = 0 ; k < s->nr_cells ; k++ ) {
            space_rebuild_recycle( s , &s->cells[k] );
            s->cells[k].sorts = NULL;
            s->cells[k].nr_tasks = 0;
            s->cells[k].nr_density = 0;
278
            s->cells[k].nr_force = 0;
279
280
            s->cells[k].density = NULL;
            s->cells[k].force = NULL;
Pedro Gonnet's avatar
Pedro Gonnet committed
281
            s->cells[k].dx_max = 0.0f;
282
            s->cells[k].sorted = 0;
Pedro Gonnet's avatar
Pedro Gonnet committed
283
            s->cells[k].count = 0;
284
            s->cells[k].gcount = 0;
285
            s->cells[k].kick1 = NULL;
Pedro Gonnet's avatar
Pedro Gonnet committed
286
            s->cells[k].kick2 = NULL;
287
            s->cells[k].super = &s->cells[k];
Pedro Gonnet's avatar
Pedro Gonnet committed
288
289
290
291
            }
        s->maxdepth = 0;
    
        }
292
        
293
294
295
296
297
298
299
300
301
302
303
304
305
    }
    

/**
 * @brief Re-build the cells as well as the tasks.
 *
 * @param s The #space in which to update the cells.
 * @param cell_max Maximal cell size.
 *
 */
 
void space_rebuild ( struct space *s , double cell_max ) {

306
    int j, k, cdim[3], nr_parts = s->nr_parts, nr_gparts = s->nr_gparts;
307
    struct cell *restrict c, *restrict cells;
308
309
    struct part *restrict finger, *restrict p, *parts = s->parts;
    struct xpart *xfinger, *xparts = s->xparts;
310
    struct gpart *gp, *gparts = s->gparts, *gfinger;
311
312
313
314
315
316
317
    int *ind;
    double ih[3], dim[3];
    // ticks tic;
    
    /* Be verbose about this. */
    // message( "re)building space..." ); fflush(stdout);
    
318
319
    /* Re-grid if necessary, or just re-set the cell data. */
    space_regrid( s , cell_max );
320
    cells = s->cells;
321
        
322
    /* Run through the particles and get their cell index. */
323
    // tic = getticks();
324
325
    const int ind_size = s->size_parts;
    if ( ( ind = (int *)malloc( sizeof(int) * ind_size ) ) == NULL )
326
        error( "Failed to allocate temporary particle indices." );
Pedro Gonnet's avatar
Pedro Gonnet committed
327
328
    ih[0] = s->ih[0]; ih[1] = s->ih[1]; ih[2] = s->ih[2];
    dim[0] = s->dim[0]; dim[1] = s->dim[1]; dim[2] = s->dim[2];
Pedro Gonnet's avatar
Pedro Gonnet committed
329
    cdim[0] = s->cdim[0]; cdim[1] = s->cdim[1]; cdim[2] = s->cdim[2];
330
    #pragma omp parallel for private(p,j)
331
332
    for ( k = 0 ; k < nr_parts ; k++ )  {
        p = &parts[k];
333
334
        for ( j = 0 ; j < 3 ; j++ )
            if ( p->x[j] < 0.0 )
Pedro Gonnet's avatar
Pedro Gonnet committed
335
336
337
338
                p->x[j] += dim[j];
            else if ( p->x[j] >= dim[j] )
                p->x[j] -= dim[j];
        ind[k] = cell_getid( cdim , p->x[0]*ih[0] , p->x[1]*ih[1] , p->x[2]*ih[2] );
339
        atomic_inc( &cells[ ind[k] ].count );
340
        }
341
342
343
344
    // message( "getting particle indices took %.3f ms." , (double)(getticks() - tic) / CPU_TPS * 1000 );


    #ifdef WITH_MPI
345
        /* Move non-local parts to the end of the list. */
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
        int nodeID = s->e->nodeID;
        for ( k = 0 ; k < nr_parts ; k++ )
            if ( cells[ ind[k] ].nodeID != nodeID ) {
                cells[ ind[k] ].count -= 1;
                nr_parts -= 1;
                struct part tp = parts[k];
                parts[k] = parts[ nr_parts ];
                parts[ nr_parts ] = tp;
                struct xpart txp = xparts[k];
                xparts[k] = xparts[ nr_parts ];
                xparts[ nr_parts ] = txp;
                int t = ind[k];
                ind[k] = ind[ nr_parts ];
                ind[ nr_parts ] = t;
                }
361
362
363
                
        /* Exchange the strays, note that this potentially re-allocates
           the parts arrays. */
364
        s->nr_parts = nr_parts + engine_exchange_strays( s->e , nr_parts , &ind[nr_parts] , s->nr_parts - nr_parts );
365
366
        parts = s->parts;
        xparts = s->xparts;
367
368
369
370
371
372
373
374
375
376
377
        
        /* Re-allocate the index array if needed.. */
        if (s->nr_parts > ind_size) {
          int *ind_new;
          if ( ( ind_new = (int *)malloc( sizeof(int) * s->nr_parts ) ) == NULL )
              error( "Failed to allocate temporary particle indices." );
          memcpy(ind_new, ind, sizeof(int) * nr_parts);
          free(ind); ind = ind_new;
        }
        
        /* Assign each particle to its cell. */
378
379
380
381
        for ( k = nr_parts ; k < s->nr_parts ; k++ ) {
            p = &parts[k];
            ind[k] = cell_getid( cdim , p->x[0]*ih[0] , p->x[1]*ih[1] , p->x[2]*ih[2] );
            cells[ ind[k] ].count += 1;
382
383
            /* if ( cells[ ind[k] ].nodeID != nodeID )
                error( "Received part that does not belong to me (nodeID=%i)." , cells[ ind[k] ].nodeID ); */
384
385
386
387
            }
        nr_parts = s->nr_parts;
    #endif
    
388
389

    /* Sort the parts according to their cells. */
390
    // tic = getticks();
391
392
    parts_sort( parts , xparts , ind , nr_parts , 0 , s->nr_cells-1 );
    // message( "parts_sort took %.3f ms." , (double)(getticks() - tic) / CPU_TPS * 1000 );
393
    
394
395
396
397
398
    /* Re-link the gparts. */
    for ( k = 0 ; k < nr_parts ; k++ )
        if ( parts[k].gpart != NULL )
            parts[k].gpart->part = &parts[k];
    
Pedro Gonnet's avatar
Pedro Gonnet committed
399
400
401
402
403
404
405
406
407
    /* Verify sort. */
    /* for ( k = 1 ; k < nr_parts ; k++ ) {
        if ( ind[k-1] > ind[k] ) {
            error( "Sort failed!" );
            }
        else if ( ind[k] != cell_getid( cdim , parts[k].x[0]*ih[0] , parts[k].x[1]*ih[1] , parts[k].x[2]*ih[2] ) )
            error( "Incorrect indices!" );
        } */
    
408
409
    /* We no longer need the indices as of here. */
    free( ind );    
410

411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433


    /* Run through the gravity particles and get their cell index. */
    // tic = getticks();
    if ( ( ind = (int *)malloc( sizeof(int) * s->size_gparts ) ) == NULL )
        error( "Failed to allocate temporary particle indices." );
    #pragma omp parallel for private(gp,j)
    for ( k = 0 ; k < nr_gparts ; k++ )  {
        gp = &gparts[k];
        for ( j = 0 ; j < 3 ; j++ )
            if ( gp->x[j] < 0.0 )
                gp->x[j] += dim[j];
            else if ( gp->x[j] >= dim[j] )
                gp->x[j] -= dim[j];
        ind[k] = cell_getid( cdim , gp->x[0]*ih[0] , gp->x[1]*ih[1] , gp->x[2]*ih[2] );
        atomic_inc( &cells[ ind[k] ].gcount );
        }
    // message( "getting particle indices took %.3f ms." , (double)(getticks() - tic) / CPU_TPS * 1000 );

    /* TODO: Here we should exchange the gparts as well! */

    /* Sort the parts according to their cells. */
    // tic = getticks();
434
    gparts_sort( gparts ,ind , nr_gparts , 0 , s->nr_cells-1 );
435
436
437
438
439
440
441
442
443
444
445
446
    // message( "gparts_sort took %.3f ms." , (double)(getticks() - tic) / CPU_TPS * 1000 );
    
    /* Re-link the parts. */
    for ( k = 0 ; k < nr_gparts ; k++ )
        if ( gparts[k].id > 0 )
            gparts[k].part->gpart = &gparts[k];

    /* We no longer need the indices as of here. */
    free( ind );    



447
    /* Hook the cells up to the parts. */
448
    // tic = getticks();
449
450
    finger = parts;
    xfinger = xparts;
451
    gfinger = gparts;
452
    for ( k = 0 ; k < s->nr_cells ; k++ ) {
453
        c = &cells[ k ];
454
        c->parts = finger;
455
        c->xparts = xfinger;
456
        c->gparts = gfinger;
457
        finger = &finger[ c->count ];
458
        xfinger = &xfinger[ c->count ];
459
        gfinger = &gfinger[ c->gcount ];
460
        }
461
    // message( "hooking up cells took %.3f ms." , (double)(getticks() - tic) / CPU_TPS * 1000 );
462
        
463
464
    /* At this point, we have the upper-level cells, old or new. Now make
       sure that the parts in each cell are ok. */
465
    // tic = getticks();
Pedro Gonnet's avatar
Pedro Gonnet committed
466
    k = 0;
467
    #pragma omp parallel shared(s,k)
Pedro Gonnet's avatar
Pedro Gonnet committed
468
    {
469
470
471
472
        if ( omp_get_thread_num() < 8 )
            while ( 1 ) {
                int myk = atomic_inc( &k );
                if ( myk < s->nr_cells )
473
                    space_split( s , &cells[myk] );
474
475
476
                else
                    break;
                }
Pedro Gonnet's avatar
Pedro Gonnet committed
477
        }
478
479
    // message( "space_split took %.3f ms." , (double)(getticks() - tic) / CPU_TPS * 1000 );
    
480
481
482
    }


483
/**
484
 * @brief Sort the particles and condensed particles according to the given indices.
485
486
 *
 * @param parts The list of #part
487
 * @param xparts The list of reduced particles
488
489
490
491
492
493
 * @param ind The indices with respect to which the parts are sorted.
 * @param N The number of parts
 * @param min Lowest index.
 * @param max highest index.
 */
 
494
void parts_sort ( struct part *parts , struct xpart *xparts , int *ind , int N , int min , int max ) {
495

496
    struct qstack {
497
        volatile int i, j, min, max;
498
        volatile int ready;
499
500
        };
    struct qstack *qstack;
501
    int qstack_size = 2*(max-min) + 10;
502
    volatile unsigned int first, last, waiting;
503
504
505
    
    int pivot;
    int i, ii, j, jj, temp_i, qid;
506
    struct part temp_p;
507
    struct xpart temp_xp;
508
509
510
511

    /* for ( int k = 0 ; k < N ; k++ )
        if ( ind[k] > max || ind[k] < min )
	    error( "ind[%i]=%i is not in [%i,%i]." , k , ind[k] , min , max ); */
512
    
513
514
515
516
    /* Allocate the stack. */
    if ( ( qstack = malloc( sizeof(struct qstack) * qstack_size ) ) == NULL )
        error( "Failed to allocate qstack." );
    
517
518
519
520
521
    /* Init the interval stack. */
    qstack[0].i = 0;
    qstack[0].j = N-1;
    qstack[0].min = min;
    qstack[0].max = max;
522
    qstack[0].ready = 1;
523
    for ( i = 1 ; i < qstack_size ; i++ )
524
525
        qstack[i].ready = 0;
    first = 0; last = 1; waiting = 1;
526
527
    
    /* Parallel bit. */
528
    #pragma omp parallel default(shared) shared(N,first,last,waiting,qstack,parts,xparts,ind,qstack_size,stderr,engine_rank) private(pivot,i,ii,j,jj,min,max,temp_i,qid,temp_xp,temp_p)
529
    {
530
    
531
        /* Main loop. */
532
        if ( omp_get_thread_num() < 8 )
533
534
        while ( waiting > 0 ) {
        
535
            /* Grab an interval off the queue. */
536
            qid = atomic_inc( &first ) % qstack_size;
537
            
538
            /* Wait for the interval to be ready. */
539
            while ( waiting > 0 && atomic_cas( &qstack[qid].ready , 1 , 1 ) != 1 );
540
541
542
543
544
            
            /* Broke loop for all the wrong reasons? */
            if ( waiting == 0 )
                break;
        
545
546
547
548
549
            /* Get the stack entry. */
            i = qstack[qid].i;
            j = qstack[qid].j;
            min = qstack[qid].min;
            max = qstack[qid].max;
550
            qstack[qid].ready = 0;
551
            // message( "thread %i got interval [%i,%i] with values in [%i,%i]." , omp_get_thread_num() , i , j , min , max );
552
            
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
            /* Loop over sub-intervals. */
            while ( 1 ) {
            
                /* Bring beer. */
                pivot = (min + max) / 2;
                
                /* One pass of QuickSort's partitioning. */
                ii = i; jj = j;
                while ( ii < jj ) {
                    while ( ii <= j && ind[ii] <= pivot )
                        ii++;
                    while ( jj >= i && ind[jj] > pivot )
                        jj--;
                    if ( ii < jj ) {
                        temp_i = ind[ii]; ind[ii] = ind[jj]; ind[jj] = temp_i;
                        temp_p = parts[ii]; parts[ii] = parts[jj]; parts[jj] = temp_p;
569
                        temp_xp = xparts[ii]; xparts[ii] = xparts[jj]; xparts[jj] = temp_xp;
570
                        }
571
                    }
572

573
574
575
                /* Verify sort. */
                /* for ( int k = i ; k <= jj ; k++ )
                    if ( ind[k] > pivot ) {
576
                        message( "sorting failed at k=%i, ind[k]=%i, pivot=%i, i=%i, j=%i, N=%i." , k , ind[k] , pivot , i , j , N );
577
578
579
580
                        error( "Partition failed (<=pivot)." );
                        }
                for ( int k = jj+1 ; k <= j ; k++ )
                    if ( ind[k] <= pivot ) {
581
                        message( "sorting failed at k=%i, ind[k]=%i, pivot=%i, i=%i, j=%i, N=%i." , k , ind[k] , pivot , i , j , N );
582
583
584
585
586
587
588
589
                        error( "Partition failed (>pivot)." );
                        } */
                        
                /* Split-off largest interval. */
                if ( jj - i > j - jj+1 ) {

                    /* Recurse on the left? */
                    if ( jj > i  && pivot > min ) {
590
                        qid = atomic_inc( &last ) % qstack_size;
591
                        while ( atomic_cas( &qstack[qid].ready , 0 , 0 ) != 0 );
592
593
594
595
596
                        qstack[qid].i = i;
                        qstack[qid].j = jj;
                        qstack[qid].min = min;
                        qstack[qid].max = pivot;
                        qstack[qid].ready = 1;
597
598
                        if ( atomic_inc( &waiting ) >= qstack_size )
                            error( "Qstack overflow." );
599
600
601
602
603
604
605
606
607
608
                        }

                    /* Recurse on the right? */
                    if ( jj+1 < j && pivot+1 < max ) {
                        i = jj+1;
                        min = pivot+1;
                        }
                    else
                        break;
                        
609
                    }
610
611
612
613
614
                    
                else {
                
                    /* Recurse on the right? */
                    if ( jj+1 < j && pivot+1 < max ) {
615
                        qid = atomic_inc( &last ) % qstack_size;
616
                        while ( atomic_cas( &qstack[qid].ready , 0 , 0 ) != 0 );
617
618
619
620
621
                        qstack[qid].i = jj+1;
                        qstack[qid].j = j;
                        qstack[qid].min = pivot+1;
                        qstack[qid].max = max;
                        qstack[qid].ready = 1;
622
623
                        if ( atomic_inc( &waiting ) >= qstack_size )
                            error( "Qstack overflow." );
624
625
626
627
628
629
630
631
632
                        }
                        
                    /* Recurse on the left? */
                    if ( jj > i  && pivot > min ) {
                        j = jj;
                        max = pivot;
                        }
                    else
                        break;
633

634
635
636
                    }
                    
                } /* loop over sub-intervals. */
637
    
638
            atomic_dec( &waiting );
639

640
641
642
            } /* main loop. */
    
        } /* parallel bit. */
643
    
644
    /* Verify sort. */
645
    /* for ( i = 1 ; i < N ; i++ )
646
        if ( ind[i-1] > ind[i] )
647
            error( "Sorting failed (ind[%i]=%i,ind[%i]=%i)." , i-1 , ind[i-1] , i , ind[i] ); */
648
649
650
            
    /* Clean up. */
    free( qstack );
651

652
653
654
    }


655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
void gparts_sort ( struct gpart *gparts , int *ind , int N , int min , int max ) {

    struct qstack {
        volatile int i, j, min, max;
        volatile int ready;
        };
    struct qstack *qstack;
    int qstack_size = 2*(max-min) + 10;
    volatile unsigned int first, last, waiting;
    
    int pivot;
    int i, ii, j, jj, temp_i, qid;
    struct gpart temp_p;

    /* for ( int k = 0 ; k < N ; k++ )
        if ( ind[k] > max || ind[k] < min )
	    error( "ind[%i]=%i is not in [%i,%i]." , k , ind[k] , min , max ); */
    
    /* Allocate the stack. */
    if ( ( qstack = malloc( sizeof(struct qstack) * qstack_size ) ) == NULL )
        error( "Failed to allocate qstack." );
    
    /* Init the interval stack. */
    qstack[0].i = 0;
    qstack[0].j = N-1;
    qstack[0].min = min;
    qstack[0].max = max;
    qstack[0].ready = 1;
    for ( i = 1 ; i < qstack_size ; i++ )
        qstack[i].ready = 0;
    first = 0; last = 1; waiting = 1;
    
    /* Parallel bit. */
688
    #pragma omp parallel default(shared) shared(N,first,last,waiting,qstack,gparts,ind,qstack_size,stderr,engine_rank) private(pivot,i,ii,j,jj,min,max,temp_i,qid,temp_p)
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
    {
    
        /* Main loop. */
        if ( omp_get_thread_num() < 8 )
        while ( waiting > 0 ) {
        
            /* Grab an interval off the queue. */
            qid = atomic_inc( &first ) % qstack_size;
            
            /* Wait for the interval to be ready. */
            while ( waiting > 0 && atomic_cas( &qstack[qid].ready , 1 , 1 ) != 1 );
            
            /* Broke loop for all the wrong reasons? */
            if ( waiting == 0 )
                break;
        
            /* Get the stack entry. */
            i = qstack[qid].i;
            j = qstack[qid].j;
            min = qstack[qid].min;
            max = qstack[qid].max;
            qstack[qid].ready = 0;
            // message( "thread %i got interval [%i,%i] with values in [%i,%i]." , omp_get_thread_num() , i , j , min , max );
            
            /* Loop over sub-intervals. */
            while ( 1 ) {
            
                /* Bring beer. */
                pivot = (min + max) / 2;
                
                /* One pass of QuickSort's partitioning. */
                ii = i; jj = j;
                while ( ii < jj ) {
                    while ( ii <= j && ind[ii] <= pivot )
                        ii++;
                    while ( jj >= i && ind[jj] > pivot )
                        jj--;
                    if ( ii < jj ) {
                        temp_i = ind[ii]; ind[ii] = ind[jj]; ind[jj] = temp_i;
                        temp_p = gparts[ii]; gparts[ii] = gparts[jj]; gparts[jj] = temp_p;
                        }
                    }

                /* Verify sort. */
                /* for ( int k = i ; k <= jj ; k++ )
                    if ( ind[k] > pivot ) {
                        message( "sorting failed at k=%i, ind[k]=%i, pivot=%i, i=%i, j=%i, N=%i." , k , ind[k] , pivot , i , j , N );
                        error( "Partition failed (<=pivot)." );
                        }
                for ( int k = jj+1 ; k <= j ; k++ )
                    if ( ind[k] <= pivot ) {
                        message( "sorting failed at k=%i, ind[k]=%i, pivot=%i, i=%i, j=%i, N=%i." , k , ind[k] , pivot , i , j , N );
                        error( "Partition failed (>pivot)." );
                        } */
                        
                /* Split-off largest interval. */
                if ( jj - i > j - jj+1 ) {

                    /* Recurse on the left? */
                    if ( jj > i  && pivot > min ) {
                        qid = atomic_inc( &last ) % qstack_size;
                        while ( atomic_cas( &qstack[qid].ready , 0 , 0 ) != 0 );
                        qstack[qid].i = i;
                        qstack[qid].j = jj;
                        qstack[qid].min = min;
                        qstack[qid].max = pivot;
                        qstack[qid].ready = 1;
                        if ( atomic_inc( &waiting ) >= qstack_size )
                            error( "Qstack overflow." );
                        }

                    /* Recurse on the right? */
                    if ( jj+1 < j && pivot+1 < max ) {
                        i = jj+1;
                        min = pivot+1;
                        }
                    else
                        break;
                        
                    }
                    
                else {
                
                    /* Recurse on the right? */
                    if ( jj+1 < j && pivot+1 < max ) {
                        qid = atomic_inc( &last ) % qstack_size;
                        while ( atomic_cas( &qstack[qid].ready , 0 , 0 ) != 0 );
                        qstack[qid].i = jj+1;
                        qstack[qid].j = j;
                        qstack[qid].min = pivot+1;
                        qstack[qid].max = max;
                        qstack[qid].ready = 1;
                        if ( atomic_inc( &waiting ) >= qstack_size )
                            error( "Qstack overflow." );
                        }
                        
                    /* Recurse on the left? */
                    if ( jj > i  && pivot > min ) {
                        j = jj;
                        max = pivot;
                        }
                    else
                        break;

                    }
                    
                } /* loop over sub-intervals. */
    
            atomic_dec( &waiting );

            } /* main loop. */
    
        } /* parallel bit. */
    
    /* Verify sort. */
    /* for ( i = 1 ; i < N ; i++ )
        if ( ind[i-1] > ind[i] )
            error( "Sorting failed (ind[%i]=%i,ind[%i]=%i)." , i-1 , ind[i-1] , i , ind[i] ); */
            
    /* Clean up. */
    free( qstack );

    }


Pedro Gonnet's avatar
Pedro Gonnet committed
814
/**
815
 * @brief Mapping function to free the sorted indices buffers.
Pedro Gonnet's avatar
Pedro Gonnet committed
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
 */

void space_map_clearsort ( struct cell *c , void *data ) {

    if ( c->sort != NULL ) {
        free( c->sort );
        c->sort = NULL;
        }

    }


/**
 * @brief Map a function to all particles in a aspace.
 *
 * @param s The #space we are working in.
832
833
 * @param fun Function pointer to apply on the cells.
 * @param data Data passed to the function fun.
Pedro Gonnet's avatar
Pedro Gonnet committed
834
835
836
837
 */
 
void space_map_parts ( struct space *s , void (*fun)( struct part *p , struct cell *c , void *data ) , void *data ) {

Pedro Gonnet's avatar
bug.    
Pedro Gonnet committed
838
    int cid = 0;
Pedro Gonnet's avatar
Pedro Gonnet committed
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857

    void rec_map ( struct cell *c ) {
    
        int k;
        
        /* No progeny? */
        if ( !c->split )
            for ( k = 0 ; k < c->count ; k++ )
                fun( &c->parts[k] , c , data );
                
        /* Otherwise, recurse. */
        else
            for ( k = 0 ; k < 8 ; k++ )
                if ( c->progeny[k] != NULL )
                    rec_map( c->progeny[k] );
                
        }
        
    /* Call the recursive function on all higher-level cells. */
858
859
860
861
862
863
864
    #pragma omp parallel shared(cid)
    {
        int mycid;
        while ( 1 ) {
            #pragma omp critical
            mycid = cid++;
            if ( mycid < s->nr_cells )
Pedro Gonnet's avatar
bug.    
Pedro Gonnet committed
865
                rec_map( &s->cells[mycid] );
866
867
868
869
            else
                break;
            }
        }
Pedro Gonnet's avatar
Pedro Gonnet committed
870
871
872
873
874
875
876
877

    }


/**
 * @brief Map a function to all particles in a aspace.
 *
 * @param s The #space we are working in.
878
 * @param full Map to all cells, including cells with sub-cells.
879
880
 * @param fun Function pointer to apply on the cells.
 * @param data Data passed to the function fun.
Pedro Gonnet's avatar
Pedro Gonnet committed
881
882
 */
 
883
884
void space_map_cells_post ( struct space *s , int full , void (*fun)( struct cell *c , void *data ) , void *data ) {

Pedro Gonnet's avatar
bug.    
Pedro Gonnet committed
885
    int cid = 0;
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903

    void rec_map ( struct cell *c ) {
    
        int k;
        
        /* Recurse. */
        if ( c->split )
            for ( k = 0 ; k < 8 ; k++ )
                if ( c->progeny[k] != NULL )
                    rec_map( c->progeny[k] );
                
        /* No progeny? */
        if ( full || !c->split )
            fun( c , data );
                
        }
        
    /* Call the recursive function on all higher-level cells. */
904
    // #pragma omp parallel shared(s,cid)
905
906
907
    {
        int mycid;
        while ( 1 ) {
908
            // #pragma omp critical
909
910
            mycid = cid++;
            if ( mycid < s->nr_cells )
Pedro Gonnet's avatar
bug.    
Pedro Gonnet committed
911
                rec_map( &s->cells[mycid] );
912
913
914
915
            else
                break;
            }
        }
916
917
918
919
920

    }


void space_map_cells_pre ( struct space *s , int full , void (*fun)( struct cell *c , void *data ) , void *data ) {
Pedro Gonnet's avatar
Pedro Gonnet committed
921

Pedro Gonnet's avatar
bug.    
Pedro Gonnet committed
922
    int cid = 0;
Pedro Gonnet's avatar
Pedro Gonnet committed
923
924
925
926
927
928

    void rec_map ( struct cell *c ) {
    
        int k;
        
        /* No progeny? */
929
        if ( full || !c->split )
Pedro Gonnet's avatar
Pedro Gonnet committed
930
931
            fun( c , data );
                
932
933
        /* Recurse. */
        if ( c->split )
Pedro Gonnet's avatar
Pedro Gonnet committed
934
935
936
937
938
939
940
            for ( k = 0 ; k < 8 ; k++ )
                if ( c->progeny[k] != NULL )
                    rec_map( c->progeny[k] );
                
        }
        
    /* Call the recursive function on all higher-level cells. */
941
    // #pragma omp parallel shared(s,cid)
Pedro Gonnet's avatar
Pedro Gonnet committed
942
943
944
    {
        int mycid;
        while ( 1 ) {
945
            // #pragma omp critical
Pedro Gonnet's avatar
Pedro Gonnet committed
946
947
            mycid = cid++;
            if ( mycid < s->nr_cells )
Pedro Gonnet's avatar
bug.    
Pedro Gonnet committed
948
                rec_map( &s->cells[mycid] );
Pedro Gonnet's avatar
Pedro Gonnet committed
949
950
951
952
            else
                break;
            }
        }
Pedro Gonnet's avatar
Pedro Gonnet committed
953
954
955
956
957
958
959
960
961
962
963
964
965

    }


/**
 * @brief Split cells that contain too many particles.
 *
 * @param s The #space we are working in.
 * @param c The #cell under consideration.
 */
 
void space_split ( struct space *s , struct cell *c ) {

966
    int k, count = c->count, gcount = c->gcount, maxdepth = 0;
967
    float h, h_max = 0.0f, dt, dt_min = FLT_MAX, dt_max = FLT_MIN;
968
969
    struct cell *temp;
    struct part *p, *parts = c->parts;
970
    struct xpart *xp, *xparts = c->xparts;
971
972
973
974
975
976
    
    /* Check the depth. */
    if ( c->depth > s->maxdepth )
        s->maxdepth = c->depth;
    
    /* Split or let it be? */
977
    if ( count > space_splitsize || gcount > space_splitsize ) {
978
979
980
981
982
983
984
985
    
        /* No longer just a leaf. */
        c->split = 1;
        
        /* Create the cell's progeny. */
        for ( k = 0 ; k < 8 ; k++ ) {
            temp = space_getcell( s );
            temp->count = 0;
986
            temp->gcount = 0;
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
            temp->loc[0] = c->loc[0];
            temp->loc[1] = c->loc[1];
            temp->loc[2] = c->loc[2];
            temp->h[0] = c->h[0]/2;
            temp->h[1] = c->h[1]/2;
            temp->h[2] = c->h[2]/2;
            temp->dmin = c->dmin/2;
            if ( k & 4 )
                temp->loc[0] += temp->h[0];
            if ( k & 2 )
                temp->loc[1] += temp->h[1];
            if ( k & 1 )
                temp->loc[2] += temp->h[2];
            temp->depth = c->depth + 1;