space.c 38 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
/*******************************************************************************
* This file is part of SWIFT.
* Copyright (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/>.
*
******************************************************************************/
Pedro Gonnet's avatar
Pedro Gonnet committed
19
20
21
22
23
24
25
26

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

/* Some standard headers. */
#include <float.h>
#include <limits.h>
#include <math.h>
27
#include <string.h>
28
#include <stdlib.h>
Pedro Gonnet's avatar
Pedro Gonnet committed
29

30
31
/* MPI headers. */
#ifdef WITH_MPI
32
#include <mpi.h>
33
34
#endif

35
36
37
/* This object's header. */
#include "space.h"

Pedro Gonnet's avatar
Pedro Gonnet committed
38
/* Local headers. */
39
#include "atomic.h"
40
#include "engine.h"
41
#include "error.h"
42
43
#include "kernel.h"
#include "lock.h"
44
#include "minmax.h"
45
#include "runner.h"
Pedro Gonnet's avatar
Pedro Gonnet committed
46

47
48
49
/* Shared sort structure. */
struct parallel_sort space_sort_struct;

Pedro Gonnet's avatar
Pedro Gonnet committed
50
51
/* Split size. */
int space_splitsize = space_splitsize_default;
52
int space_subsize = space_subsize_default;
53
int space_maxsize = space_maxsize_default;
Pedro Gonnet's avatar
Pedro Gonnet committed
54
55
56

/* Map shift vector to sortlist. */
const int sortlistID[27] = {
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
    /* ( -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};

85
86
87
88
89
90
91
92
93
94
95
96
/**
 * @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.
 */

97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
int space_getsid(struct space *s, struct cell **ci, struct cell **cj,
                 double *shift) {

  int k, sid = 0, periodic = s->periodic;
  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];
    if (periodic && dx[k] < -s->dim[k] / 2)
      shift[k] = s->dim[k];
    else if (periodic && dx[k] > s->dim[k] / 2)
      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;
}
132

133
/**
134
 * @brief Recursively dismantle a cell tree.
135
136
 *
 */
137
138
139
140
141
142
143
144
145
146
147
148
149
150

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;
      }
}

151
/**
152
 * @brief Re-build the cell grid.
153
 *
154
155
 * @param s The #space.
 * @param cell_max Maximum cell edge length.
156
 * @param verbose Print messages to stdout or not.
157
 */
158

159
void space_regrid(struct space *s, double cell_max, int verbose) {
160
161
162
163

  float h_max = s->cell_min / kernel_gamma / space_stretch, dmin;
  int i, j, k, cdim[3], nr_parts = s->nr_parts;
  struct cell *restrict c;
164
  ticks tic = getticks();
165
166
167
168
169
170

  /* Run through the parts and get the current h_max. */
  // tic = getticks();
  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;
171
    }
172
173
174
175
176
177
178
179
180
181
182
183
184
185
  } 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;
  }

/* 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)
186
      error("Failed to aggregate the rebuild flag across nodes.");
187
188
189
    h_max = buff;
  }
#endif
190
  if (verbose) message("h_max is %.3e (cell_max=%.3e).", h_max, cell_max);
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257

  /* Get the new putative cell dimensions. */
  for (k = 0; k < 3; k++)
    cdim[k] =
        floor(s->dim[k] / fmax(h_max * kernel_gamma * space_stretch, cell_max));

  /* 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.");

/* 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

  /* Do we need to re-build the upper-level cells? */
  // tic = getticks();
  if (s->cells == NULL || cdim[0] < s->cdim[0] || cdim[1] < s->cdim[1] ||
      cdim[2] < s->cdim[2]) {

    /* Free the old cells, if they were allocated. */
    if (s->cells != NULL) {
      for (k = 0; k < s->nr_cells; k++) {
        space_rebuild_recycle(s, &s->cells[k]);
        if (s->cells[k].sort != NULL) free(s->cells[k].sort);
      }
      free(s->cells);
      s->maxdepth = 0;
    }

    /* Set the new cell dimensions only if smaller. */
    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];
    }
    dmin = fminf(s->h[0], fminf(s->h[1], s->h[2]));

    /* Allocate the highest level of cells. */
    s->tot_cells = s->nr_cells = cdim[0] * cdim[1] * cdim[2];
    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];
          c->dmin = dmin;
          c->depth = 0;
          c->count = 0;
          c->gcount = 0;
          c->super = c;
          lock_init(&c->lock);
Pedro Gonnet's avatar
Pedro Gonnet committed
258
        }
259
260

    /* Be verbose about the change. */
261
262
263
    if (verbose)
      message("set cell dimensions to [ %i %i %i ].", cdim[0], cdim[1],
              cdim[2]);
264
265
266
    fflush(stdout);

  } /* re-build upper-level cells? */
267
268
  // message( "rebuilding upper-level cells took %.3f %s." ,
  // clocks_from_ticks(double)(getticks() - tic), clocks_getunit());
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285

  /* 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;
      s->cells[k].nr_force = 0;
      s->cells[k].density = NULL;
      s->cells[k].force = NULL;
      s->cells[k].dx_max = 0.0f;
      s->cells[k].sorted = 0;
      s->cells[k].count = 0;
      s->cells[k].gcount = 0;
Matthieu Schaller's avatar
Matthieu Schaller committed
286
      s->cells[k].init = NULL;
Matthieu Schaller's avatar
Matthieu Schaller committed
287
      s->cells[k].ghost = NULL;
Matthieu Schaller's avatar
Matthieu Schaller committed
288
289
      s->cells[k].drift = NULL;
      s->cells[k].kick = NULL;
290
      s->cells[k].super = &s->cells[k];
291
    }
292
293
    s->maxdepth = 0;
  }
294
295
296
297

  if (verbose)
    message("took %.3f %s.", clocks_from_ticks(getticks() - tic),
            clocks_getunit());
298
}
299
300
301
302
303
304

/**
 * @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.
305
 * @param verbose Print messages to stdout or not
306
307
 *
 */
308

309
void space_rebuild(struct space *s, double cell_max, int verbose) {
310

311
  ticks tic = getticks();
312
313
314
315
316

  /* Be verbose about this. */
  // message( "re)building space..." ); fflush(stdout);

  /* Re-grid if necessary, or just re-set the cell data. */
317
  space_regrid(s, cell_max, verbose);
318

319
320
321
322
323
324
  int nr_parts = s->nr_parts;
  int nr_gparts = s->nr_gparts;
  struct cell *restrict cells = s->cells;

  double ih[3], dim[3];
  int cdim[3];
325
326
327
328
329
330
331
332
333
  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];
  cdim[0] = s->cdim[0];
  cdim[1] = s->cdim[1];
  cdim[2] = s->cdim[2];
334
335
336
337
338
339
340
341
342
343

  /* Run through the particles and get their cell index. */
  // tic = getticks();
  const size_t ind_size = s->size_parts;
  size_t *ind;
  if ((ind = (size_t *)malloc(sizeof(size_t) * ind_size)) == NULL)
    error("Failed to allocate temporary particle indices.");
  for (int k = 0; k < nr_parts; k++) {
    struct part *restrict p = &s->parts[k];
    for (int j = 0; j < 3; j++)
344
345
346
347
      if (p->x[j] < 0.0)
        p->x[j] += dim[j];
      else if (p->x[j] >= dim[j])
        p->x[j] -= dim[j];
348
    ind[k] =
349
        cell_getid(cdim, p->x[0] * ih[0], p->x[1] * ih[1], p->x[2] * ih[2]);
350
    cells[ind[k]].count++;
351
  }
Matthieu Schaller's avatar
Matthieu Schaller committed
352
353
// message( "getting particle indices took %.3f %s." ,
// clocks_from_ticks(getticks() - tic), clocks_getunit()):
354
355
356

#ifdef WITH_MPI
  /* Move non-local parts to the end of the list. */
357
358
  const int nodeID = s->e->nodeID;
  for (int k = 0; k < nr_parts; k++)
359
360
361
    if (cells[ind[k]].nodeID != nodeID) {
      cells[ind[k]].count -= 1;
      nr_parts -= 1;
362
363
364
365
366
367
      struct part tp = s->parts[k];
      s->parts[k] = s->parts[nr_parts];
      s->parts[nr_parts] = tp;
      struct xpart txp = s->xparts[k];
      s->xparts[k] = s->xparts[nr_parts];
      s->xparts[nr_parts] = txp;
368
369
370
      int t = ind[k];
      ind[k] = ind[nr_parts];
      ind[nr_parts] = t;
371
372
    }

373
374
375
  /* Exchange the strays, note that this potentially re-allocates
     the parts arrays. */
  s->nr_parts =
376
377
      nr_parts + engine_exchange_strays(s->e, nr_parts, &ind[nr_parts],
                                        s->nr_parts - nr_parts);
378
379

  /* Re-allocate the index array if needed.. */
380
  if (s->nr_parts > ind_size) {
381
382
    size_t *ind_new;
    if ((ind_new = (size_t *)malloc(sizeof(size_t) * s->nr_parts)) == NULL)
383
      error("Failed to allocate temporary particle indices.");
384
    memcpy(ind_new, ind, sizeof(size_t) * nr_parts);
385
386
    free(ind);
    ind = ind_new;
387
388
389
  }

  /* Assign each particle to its cell. */
390
391
  for (int k = nr_parts; k < s->nr_parts; k++) {
    struct part *p = &s->parts[k];
392
    ind[k] =
393
        cell_getid(cdim, p->x[0] * ih[0], p->x[1] * ih[1], p->x[2] * ih[2]);
394
395
396
397
    cells[ind[k]].count += 1;
    /* if ( cells[ ind[k] ].nodeID != nodeID )
        error( "Received part that does not belong to me (nodeID=%i)." , cells[
       ind[k] ].nodeID ); */
398
  }
399
  nr_parts = s->nr_parts;
400
401
402
#endif

  /* Sort the parts according to their cells. */
403
  space_parts_sort(s, ind, nr_parts, 0, s->nr_cells - 1, verbose);
404
405

  /* Re-link the gparts. */
406
  part_relink_gparts(s->parts, nr_parts, 0);
407

408
  /* Verify space_sort_struct. */
409
  /* for ( k = 1 ; k < nr_parts ; k++ ) {
410
      if ( ind[k-1] > ind[k] ) {
411
412
          error( "Sort failed!" );
          }
413
      else if ( ind[k] != cell_getid( cdim , parts[k].x[0]*ih[0] ,
414
415
416
417
418
     parts[k].x[1]*ih[1] , parts[k].x[2]*ih[2] ) )
          error( "Incorrect indices!" );
      } */

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

  /* Run through the gravity particles and get their cell index. */
  // tic = getticks();
423
424
425
426
427
  const size_t gind_size = s->size_gparts;
  size_t *gind;
  if ((gind = (size_t *)malloc(sizeof(size_t) * gind_size)) == NULL)
    error("Failed to allocate temporary g-particle indices.");
  for (int k = 0; k < nr_gparts; k++) {
428
    struct gpart *gp = &s->gparts[k];
429
    for (int j = 0; j < 3; j++)
430
431
432
433
      if (gp->x[j] < 0.0)
        gp->x[j] += dim[j];
      else if (gp->x[j] >= dim[j])
        gp->x[j] -= dim[j];
434
    gind[k] =
435
        cell_getid(cdim, gp->x[0] * ih[0], gp->x[1] * ih[1], gp->x[2] * ih[2]);
436
    cells[gind[k]].gcount++;
437
  }
438
439
// message( "getting particle indices took %.3f %s." ,
// clocks_from_ticks(getticks() - tic), clocks_getunit());
440

441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
#ifdef WITH_MPI

  /* Move non-local gparts to the end of the list. */
  for (int k = 0; k < nr_gparts; k++)
    if (cells[ind[k]].nodeID != nodeID) {
      cells[ind[k]].gcount -= 1;
      nr_gparts -= 1;
      struct gpart tp = s->gparts[k];
      s->gparts[k] = s->gparts[nr_gparts];
      s->gparts[nr_gparts] = tp;
      int t = ind[k];
      ind[k] = ind[nr_gparts];
      ind[nr_gparts] = t;
    }

  /* Exchange the strays, note that this potentially re-allocates
     the parts arrays. */
Matthieu Schaller's avatar
Matthieu Schaller committed
458
  // s->nr_gparts =
459
460
  //    nr_gparts + engine_exchange_strays(s->e, nr_gparts, &ind[nr_gparts],
  //                                        s->nr_gparts - nr_gparts);
Matthieu Schaller's avatar
Matthieu Schaller committed
461
  if (nr_gparts > 0)
462
463
    error("Need to implement the exchange of strays for the gparts");

464
465
466
467
468
469
470
471
472
473
474
  /* Re-allocate the index array if needed.. */
  if (s->nr_gparts > gind_size) {
    size_t *gind_new;
    if ((gind_new = (size_t *)malloc(sizeof(size_t) * s->nr_gparts)) == NULL)
      error("Failed to allocate temporary g-particle indices.");
    memcpy(gind_new, gind, sizeof(size_t) * nr_gparts);
    free(gind);
    gind = gind_new;
  }

  /* Assign each particle to its cell. */
475
  for (int k = nr_gparts; k < s->nr_gparts; k++) {
476
477
478
479
480
481
482
483
484
485
486
    struct gpart *p = &s->gparts[k];
    gind[k] =
        cell_getid(cdim, p->x[0] * ih[0], p->x[1] * ih[1], p->x[2] * ih[2]);
    cells[gind[k]].count += 1;
    /* if ( cells[ ind[k] ].nodeID != nodeID )
        error( "Received part that does not belong to me (nodeID=%i)." , cells[
       ind[k] ].nodeID ); */
  }
  nr_gparts = s->nr_gparts;

#endif
487
488

  /* Sort the parts according to their cells. */
489
  space_gparts_sort(s->gparts, gind, nr_gparts, 0, s->nr_cells - 1);
490
491

  /* Re-link the parts. */
492
  part_relink_parts(s->gparts, nr_gparts, s->parts);
493
494

  /* We no longer need the indices as of here. */
495
  free(gind);
496
497
498

  /* Hook the cells up to the parts. */
  // tic = getticks();
499
500
501
  struct part *finger = s->parts;
  struct xpart *xfinger = s->xparts;
  struct gpart *gfinger = s->gparts;
502
503
  for (int k = 0; k < s->nr_cells; k++) {
    struct cell *restrict c = &cells[k];
504
505
506
507
508
509
510
    c->parts = finger;
    c->xparts = xfinger;
    c->gparts = gfinger;
    finger = &finger[c->count];
    xfinger = &xfinger[c->count];
    gfinger = &gfinger[c->gcount];
  }
511
  // message( "hooking up cells took %.3f %s." ,
Matthieu Schaller's avatar
Matthieu Schaller committed
512
  // clocks_from_ticks(getticks() - tic), clocks_getunit());
513
514
515

  /* At this point, we have the upper-level cells, old or new. Now make
     sure that the parts in each cell are ok. */
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
  space_split(s, cells, verbose);

  if (verbose)
    message("took %.3f %s.", clocks_from_ticks(getticks() - tic),
            clocks_getunit());
}

/**
 * @brief Split particles between cells of a hierarchy
 *
 * @param s The #space.
 * @param cells The cell hierarchy
 * @param verbose Are we talkative ?
 */
void space_split(struct space *s, struct cell *cells, int verbose) {

  ticks tic = getticks();

  for (int k = 0; k < s->nr_cells; k++)
535
536
    scheduler_addtask(&s->e->sched, task_type_split_cell, task_subtype_none, k,
                      0, &cells[k], NULL, 0);
537
  engine_launch(s->e, s->e->nr_threads, 1 << task_type_split_cell, 0);
538

539
540
541
  if (verbose)
    message("took %.3f %s.", clocks_from_ticks(getticks() - tic),
            clocks_getunit());
542
}
543

544
/**
545
546
 * @brief Sort the particles and condensed particles according to the given
 *indices.
547
 *
Matthieu Schaller's avatar
Matthieu Schaller committed
548
 * @param s The #space.
549
550
551
552
 * @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.
553
 * @param verbose Are we talkative ?
554
 */
555

556
void space_parts_sort(struct space *s, size_t *ind, size_t N, int min, int max,
557
558
559
560
561
                      int verbose) {

  ticks tic = getticks();

  /*Populate the global parallel_sort structure with the input data */
562
563
564
  space_sort_struct.parts = s->parts;
  space_sort_struct.xparts = s->xparts;
  space_sort_struct.ind = ind;
565
  space_sort_struct.stack_size = 2 * (max - min + 1) + 10 + s->e->nr_threads;
566
567
568
569
570
571
  if ((space_sort_struct.stack = malloc(sizeof(struct qstack) *
                                        space_sort_struct.stack_size)) == NULL)
    error("Failed to allocate sorting stack.");
  for (int i = 0; i < space_sort_struct.stack_size; i++)
    space_sort_struct.stack[i].ready = 0;

572
  /* Add the first interval. */
573
574
575
576
577
578
579
580
581
  space_sort_struct.stack[0].i = 0;
  space_sort_struct.stack[0].j = N - 1;
  space_sort_struct.stack[0].min = min;
  space_sort_struct.stack[0].max = max;
  space_sort_struct.stack[0].ready = 1;
  space_sort_struct.first = 0;
  space_sort_struct.last = 1;
  space_sort_struct.waiting = 1;

582
  /* Launch the sorting tasks. */
583
  engine_launch(s->e, s->e->nr_threads, (1 << task_type_psort), 0);
584
585

  /* Verify space_sort_struct. */
586
  /* for (int i = 1; i < N; i++)
587
    if (ind[i - 1] > ind[i])
588
589
      error("Sorting failed (ind[%i]=%i,ind[%i]=%i), min=%i, max=%i.", i - 1,
  ind[i - 1], i,
590
591
            ind[i], min, max);
  message("Sorting succeeded."); */
592

593
  /* Clean up. */
594
  free(space_sort_struct.stack);
595
596
597
598

  if (verbose)
    message("took %.3f %s.", clocks_from_ticks(getticks() - tic),
            clocks_getunit());
599
}
600

601
void space_do_parts_sort() {
602

603
  /* Pointers to the sorting data. */
604
  size_t *ind = space_sort_struct.ind;
605
606
  struct part *parts = space_sort_struct.parts;
  struct xpart *xparts = space_sort_struct.xparts;
607

608
  /* Main loop. */
609
  while (space_sort_struct.waiting) {
610

611
    /* Grab an interval off the queue. */
612
613
    int qid =
        atomic_inc(&space_sort_struct.first) % space_sort_struct.stack_size;
614

615
    /* Wait for the entry to be ready, or for the sorting do be done. */
616
617
    while (!space_sort_struct.stack[qid].ready)
      if (!space_sort_struct.waiting) return;
618

619
    /* Get the stack entry. */
620
621
    ptrdiff_t i = space_sort_struct.stack[qid].i;
    ptrdiff_t j = space_sort_struct.stack[qid].j;
622
623
    int min = space_sort_struct.stack[qid].min;
    int max = space_sort_struct.stack[qid].max;
624
    space_sort_struct.stack[qid].ready = 0;
625

626
627
    /* Loop over sub-intervals. */
    while (1) {
628

629
      /* Bring beer. */
630
      const int pivot = (min + max) / 2;
631
632
      /* message("Working on interval [%i,%i] with min=%i, max=%i, pivot=%i.",
              i, j, min, max, pivot); */
633
634

      /* One pass of QuickSort's partitioning. */
635
636
      ptrdiff_t ii = i;
      ptrdiff_t jj = j;
637
638
639
640
      while (ii < jj) {
        while (ii <= j && ind[ii] <= pivot) ii++;
        while (jj >= i && ind[jj] > pivot) jj--;
        if (ii < jj) {
641
          size_t temp_i = ind[ii];
642
643
          ind[ii] = ind[jj];
          ind[jj] = temp_i;
644
          struct part temp_p = parts[ii];
645
646
          parts[ii] = parts[jj];
          parts[jj] = temp_p;
647
          struct xpart temp_xp = xparts[ii];
648
649
650
651
          xparts[ii] = xparts[jj];
          xparts[jj] = temp_xp;
        }
      }
652

653
      /* Verify space_sort_struct. */
654
655
656
657
658
659
660
661
662
663
664
665
      /* 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.", k,
                  ind[k], pivot, i, j);
          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.", k,
                  ind[k], pivot, i, j);
          error("Partition failed (>pivot).");
        } */
666
667
668
669
670
671

      /* Split-off largest interval. */
      if (jj - i > j - jj + 1) {

        /* Recurse on the left? */
        if (jj > i && pivot > min) {
672
673
          qid = atomic_inc(&space_sort_struct.last) %
                space_sort_struct.stack_size;
674
675
          while (space_sort_struct.stack[qid].ready)
            ;
676
677
678
679
680
681
682
          space_sort_struct.stack[qid].i = i;
          space_sort_struct.stack[qid].j = jj;
          space_sort_struct.stack[qid].min = min;
          space_sort_struct.stack[qid].max = pivot;
          if (atomic_inc(&space_sort_struct.waiting) >=
              space_sort_struct.stack_size)
            error("Qstack overflow.");
683
          space_sort_struct.stack[qid].ready = 1;
684
        }
685

686
687
688
689
690
691
692
693
694
695
        /* Recurse on the right? */
        if (jj + 1 < j && pivot + 1 < max) {
          i = jj + 1;
          min = pivot + 1;
        } else
          break;

      } else {

        /* Recurse on the right? */
696
        if (pivot + 1 < max) {
697
698
          qid = atomic_inc(&space_sort_struct.last) %
                space_sort_struct.stack_size;
699
700
          while (space_sort_struct.stack[qid].ready)
            ;
701
702
703
704
705
706
707
          space_sort_struct.stack[qid].i = jj + 1;
          space_sort_struct.stack[qid].j = j;
          space_sort_struct.stack[qid].min = pivot + 1;
          space_sort_struct.stack[qid].max = max;
          if (atomic_inc(&space_sort_struct.waiting) >=
              space_sort_struct.stack_size)
            error("Qstack overflow.");
708
          space_sort_struct.stack[qid].ready = 1;
709
        }
710

711
712
713
714
715
716
717
        /* Recurse on the left? */
        if (jj > i && pivot > min) {
          j = jj;
          max = pivot;
        } else
          break;
      }
718

719
720
    } /* loop over sub-intervals. */

721
    atomic_dec(&space_sort_struct.waiting);
722
723

  } /* main loop. */
724
725
}

726
void space_gparts_sort(struct gpart *gparts, size_t *ind, size_t N, int min,
727
                       int max) {
728
729

  struct qstack {
730
731
    volatile size_t i, j;
    volatile int min, max;
732
733
734
735
736
737
738
    volatile int ready;
  };
  struct qstack *qstack;
  int qstack_size = 2 * (max - min) + 10;
  volatile unsigned int first, last, waiting;

  int pivot;
739
  ptrdiff_t i, ii, j, jj, temp_i;
740
  int qid;
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
  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;

762
763
  /* Main loop. */
  while (waiting > 0) {
764

765
766
    /* Grab an interval off the queue. */
    qid = (first++) % qstack_size;
767

768
769
770
771
772
773
774
775
776
    /* Get the stack entry. */
    i = qstack[qid].i;
    j = qstack[qid].j;
    min = qstack[qid].min;
    max = qstack[qid].max;
    qstack[qid].ready = 0;

    /* Loop over sub-intervals. */
    while (1) {
777

778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
      /* 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;
        }
      }
796

797
      /* Verify space_sort_struct. */
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
      /* 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 = (last++) % qstack_size;
          qstack[qid].i = i;
          qstack[qid].j = jj;
          qstack[qid].min = min;
          qstack[qid].max = pivot;
          qstack[qid].ready = 1;
          if ((waiting++) >= qstack_size) error("Qstack overflow.");
        }
824

825
826
827
828
829
830
        /* Recurse on the right? */
        if (jj + 1 < j && pivot + 1 < max) {
          i = jj + 1;
          min = pivot + 1;
        } else
          break;
831

832
833
834
      } else {

        /* Recurse on the right? */
835
        if (pivot + 1 < max) {
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
          qid = (last++) % qstack_size;
          qstack[qid].i = jj + 1;
          qstack[qid].j = j;
          qstack[qid].min = pivot + 1;
          qstack[qid].max = max;
          qstack[qid].ready = 1;
          if ((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. */

    waiting--;

  } /* main loop. */
858

859
  /* Verify space_sort_struct. */
860
861
862
863
864
865
866
867
  /* 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);
}
868

Pedro Gonnet's avatar
Pedro Gonnet committed
869
/**
870
 * @brief Mapping function to free the sorted indices buffers.
Pedro Gonnet's avatar
Pedro Gonnet committed
871
872
 */

873
void space_map_clearsort(struct cell *c, void *data) {
Pedro Gonnet's avatar
Pedro Gonnet committed
874

875
876
877
878
879
  if (c->sort != NULL) {
    free(c->sort);
    c->sort = NULL;
  }
}
Pedro Gonnet's avatar
Pedro Gonnet committed
880

881
882
883
/**
 * @brief Map a function to all particles in a cell recursively.
 *
884
 * @param c The #cell we are working in.
885
886
887
888
 * @param fun Function pointer to apply on the cells.
 * @param data Data passed to the function fun.
 */

Pedro Gonnet's avatar
Pedro Gonnet committed
889
890
891
892
static void rec_map_parts(struct cell *c,
                          void (*fun)(struct part *p, struct cell *c,
                                      void *data),
                          void *data) {
893
894
895
896
897
898

  int k;

  /* No progeny? */
  if (!c->split)
    for (k = 0; k < c->count; k++) fun(&c->parts[k], c, data);
Pedro Gonnet's avatar
Pedro Gonnet committed
899

900
901
902
903
904
905
  /* Otherwise, recurse. */
  else
    for (k = 0; k < 8; k++)
      if (c->progeny[k] != NULL) rec_map_parts(c->progeny[k], fun, data);
}

Pedro Gonnet's avatar
Pedro Gonnet committed
906
/**
907
 * @brief Map a function to all particles in a space.
Pedro Gonnet's avatar
Pedro Gonnet committed
908
909
 *
 * @param s The #space we are working in.
910
911
 * @param fun Function pointer to apply on the cells.
 * @param data Data passed to the function fun.
Pedro Gonnet's avatar
Pedro Gonnet committed
912
913
 */

914
915
916
void space_map_parts(struct space *s,
                     void (*fun)(struct part *p, struct cell *c, void *data),
                     void *data) {
Pedro Gonnet's avatar
Pedro Gonnet committed
917

918
919
  int cid = 0;

920
  /* Call the recursive function on all higher-level cells. */
Pedro Gonnet's avatar
Pedro Gonnet committed
921
922
  for (cid = 0; cid < s->nr_cells; cid++)
    rec_map_parts(&s->cells[cid], fun, data);
923
}
924

925
926
927
928
929
930
931
932
/**
 * @brief Map a function to all particles in a cell recursively.
 *
 * @param c The #cell we are working in.
 * @param fun Function pointer to apply on the cells.
 */

static void rec_map_parts_xparts(struct cell *c,
933
934
                                 void (*fun)(struct part *p, struct xpart *xp,
                                             struct cell *c)) {
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955

  int k;

  /* No progeny? */
  if (!c->split)
    for (k = 0; k < c->count; k++) fun(&c->parts[k], &c->xparts[k], c);

  /* Otherwise, recurse. */
  else
    for (k = 0; k < 8; k++)
      if (c->progeny[k] != NULL) rec_map_parts_xparts(c->progeny[k], fun);
}

/**
 * @brief Map a function to all particles (#part and #xpart) in a space.
 *
 * @param s The #space we are working in.
 * @param fun Function pointer to apply on the particles in the cells.
 */

void space_map_parts_xparts(struct space *s,
956
957
                            void (*fun)(struct part *p, struct xpart *xp,
                                        struct cell *c)) {
958
959
960
961
962
963
964
965

  int cid = 0;

  /* Call the recursive function on all higher-level cells. */
  for (cid = 0; cid < s->nr_cells; cid++)
    rec_map_parts_xparts(&s->cells[cid], fun);
}

966
967
968
/**
 * @brief Map a function to all particles in a cell recursively.
 *
969
 * @param c The #cell we are working in.
970
971
972
973
 * @param full Map to all cells, including cells with sub-cells.
 * @param fun Function pointer to apply on the cells.
 * @param data Data passed to the function fun.
 */
974

Pedro Gonnet's avatar
Pedro Gonnet committed
975
976
977
static void rec_map_cells_post(struct cell *c, int full,
                               void (*fun)(struct cell *c, void *data),
                               void *data) {
978

979
980
981
982
983
  int k;

  /* Recurse. */
  if (c->split)
    for (k = 0; k < 8; k++)
Pedro Gonnet's avatar
Pedro Gonnet committed
984
985
986
      if (c->progeny[k] != NULL)
        rec_map_cells_post(c->progeny[k], full, fun, data);

987
988
  /* No progeny? */
  if (full || !c->split) fun(c, data);
989
}
Pedro Gonnet's avatar
Pedro Gonnet committed
990
991

/**
992
 * @brief Map a function to all particles in a aspace.
Pedro Gonnet's avatar
Pedro Gonnet committed
993
994
 *
 * @param s The #space we are working in.
995
 * @param full Map to all cells, including cells with sub-cells.
996
997
 * @param fun Function pointer to apply on the cells.
 * @param data Data passed to the function fun.
Pedro Gonnet's avatar
Pedro Gonnet committed
998
 */
999

1000
void space_map_cells_post(struct space *s, int full,
Pedro Gonnet's avatar
Pedro Gonnet committed
1001
                          void (*fun)(struct cell *c, void *data), void *data) {
1002

1003
  int cid = 0;
1004

1005
  /* Call the recursive function on all higher-level cells. */
Pedro Gonnet's avatar
Pedro Gonnet committed
1006
1007
  for (cid = 0; cid < s->nr_cells; cid++)
    rec_map_cells_post(&s->cells[cid], full, fun, data);
1008
}
1009

Pedro Gonnet's avatar
Pedro Gonnet committed
1010
1011
1012
static void rec_map_cells_pre(struct cell *c, int full,
                              void (*fun)(struct cell *c, void *data),
                              void *data) {
1013

1014
  int k;
Pedro Gonnet's avatar
Pedro Gonnet committed
1015

1016
1017
  /* No progeny? */
  if (full || !c->split) fun(c, data);
Pedro Gonnet's avatar
Pedro Gonnet committed
1018

1019
1020
1021
  /* Recurse. */
  if (c->split)
    for (k = 0; k < 8; k++)
Pedro Gonnet's avatar
Pedro Gonnet committed
1022
1023
      if (c->progeny[k] != NULL)
        rec_map_cells_pre(c->progeny[k], full, fun, data);
1024
}
Pedro Gonnet's avatar
Pedro Gonnet committed
1025

1026
1027
1028
1029
1030
1031
1032
1033
/**
 * @brief Calls function fun on the cells in the space s
 *
 * @param s The #space
 * @param full If true calls the function on all cells and not just on leaves
 * @param fun The function to call.
 * @param data Additional data passed to fun() when called
 */
1034
void space_map_cells_pre(struct space *s, int full,
Pedro Gonnet's avatar
Pedro Gonnet committed
1035
                         void (*fun)(struct cell *c, void *data), void *data) {
1036

1037
  int cid = 0;
1038
1039

  /* Call the recursive function on all higher-level cells. */
Pedro Gonnet's avatar
Pedro Gonnet committed
1040
1041
  for (cid = 0; cid < s->nr_cells; cid++)
    rec_map_cells_pre(&s->cells[cid], full, fun, data);
1042
}
Pedro Gonnet's avatar
Pedro Gonnet committed
1043
1044
1045
1046
1047
1048
1049

/**
 * @brief Split cells that contain too many particles.
 *
 * @param s The #space we are working in.
 * @param c The #cell under consideration.
 */
1050

1051
void space_do_split(struct space *s, struct cell *c) {
1052

1053
1054
1055
  const int count = c->count;
  const int gcount = c->gcount;
  int maxdepth = 0;
1056
1057
  float h_max = 0.0f;
  int ti_end_min = max_nr_timesteps, ti_end_max = 0;
1058
  struct cell *temp;
1059
1060
  struct part *parts = c->parts;
  struct gpart *gparts = c->gparts;
1061
  struct xpart *xparts = c->xparts;
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072

  /* Check the depth. */
  if (c->depth > s->maxdepth) s->maxdepth = c->depth;

  /* Split or let it be? */
  if (count > space_splitsize || gcount > space_splitsize) {

    /* No longer just a leaf. */
    c->split = 1;

    /* Create the cell's progeny. */
1073
    for (int k = 0; k < 8; k++) {
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
      temp = space_getcell(s);
      temp->count = 0;
      temp->gcount = 0;
      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;
      temp->split = 0;
      temp->h_max = 0.0;
      temp->dx_max = 0.0;
      temp->nodeID = c->nodeID;
      temp->parent = c;
      c->progeny[k] = temp;
    }

    /* Split the cell data. */
1097
    cell_split(c, c->parts - s->parts);
1098
1099

    /* Remove any progeny with zero parts. */
1100
    for (int k = 0; k < 8; k++)
1101
1102
1103
1104
      if (c->progeny[k]->count == 0 && c->progeny[k]->gcount == 0) {
        space_recycle(s, c->progeny[k]);
        c->progeny[k] = NULL;
      } else {
1105
        space_do_split(s, c->progeny[k]);
1106
        h_max = fmaxf(h_max, c->progeny[k]->h_max);
1107
1108
        ti_end_min = min(ti_end_min, c->progeny[k]->ti_end_min);
        ti_end_max = max(ti_end_max, c->progeny[k]->ti_end_max);
1109
1110
1111
1112
1113
1114
        if (c->progeny[k]->maxdepth > maxdepth)
          maxdepth = c->progeny[k]->maxdepth;
      }

    /* Set the values for this cell. */
    c->h_max = h_max;
1115
1116
    c->ti_end_min = ti_end_min;
    c->ti_end_max = ti_end_max;
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
    c->maxdepth = maxdepth;

  }

  /* Otherwise, collect the data for this cell. */
  else {

    /* Clear the progeny. */
    bzero(c->progeny, sizeof(struct cell *) * 8);
    c->split = 0;
    c->maxdepth = c->depth;

    /* Get dt_min/dt_max. */
1130
    for (int k = 0; k < count; k++) {
1131
      struct part *p = &parts[k];
1132
      struct xpart *xp = &xparts[k];
1133
1134
      const float h = p->h;
      const int ti_end = p->ti_end;
1135
1136
1137
      xp->x_old[0] = p->x[0];
      xp->x_old[1] = p->x[1];
      xp->x_old[2] = p->x[2];
1138
      if (h > h_max) h_max = h;
1139
1140
      if (ti_end < ti_end_min) ti_end_min = ti_end;
      if (ti_end > ti_end_max) ti_end_max = ti_end;
1141
    }
1142
1143
1144
1145
1146
1147
    for (int k = 0; k < gcount; k++) {
      struct gpart *p = &gparts[k];
      const int ti_end = p->ti_end;
      if (ti_end < ti_end_min) ti_end_min = ti_end;
      if (ti_end > ti_end_max) ti_end_max = ti_end;
    }
1148
    c->h_max = h_max;
1149
1150
    c->ti_end_min = ti_end_min;
    c->ti_end_max = ti_end_max;
1151
  }
1152

1153
  /* Set ownership according to the start of the parts array. */
Matthieu Schaller's avatar
Matthieu Schaller committed
1154
  if (s->nr_parts > 0)
1155
1156
    c->owner =
        ((c->parts - s->parts) % s->nr_parts) * s->nr_queues / s->nr_parts;
1157
  else