space.c 35.5 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
44
#include "kernel.h"
#include "lock.h"
#include "runner.h"
Pedro Gonnet's avatar
Pedro Gonnet committed
45

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

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

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

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

96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
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;
}
131

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

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

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

158
void space_regrid(struct space *s, double cell_max, int verbose) {
159
160
161
162
163
164
165
166
167
168
169

  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;
  // ticks tic;

  /* 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;
170
    }
171
172
173
174
175
176
177
178
179
180
181
182
183
184
  } 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)
185
      error("Failed to aggregate the rebuild flag across nodes.");
186
187
188
    h_max = buff;
  }
#endif
189
  if (verbose) message("h_max is %.3e (cell_max=%.3e).", h_max, cell_max);
190
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

  /* 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
257
        }
258
259

    /* Be verbose about the change. */
260
261
262
    if (verbose)
      message("set cell dimensions to [ %i %i %i ].", cdim[0], cdim[1],
              cdim[2]);
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
    fflush(stdout);

  } /* re-build upper-level cells? */
  // message( "rebuilding upper-level cells took %.3f ms." , (double)(getticks()
  // - tic) / CPU_TPS * 1000 );

  /* 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
285
      s->cells[k].init = NULL;
Matthieu Schaller's avatar
Matthieu Schaller committed
286
      s->cells[k].ghost = NULL;
Matthieu Schaller's avatar
Matthieu Schaller committed
287
288
      s->cells[k].drift = NULL;
      s->cells[k].kick = NULL;
289
      s->cells[k].super = &s->cells[k];
290
    }
291
292
293
    s->maxdepth = 0;
  }
}
294
295
296
297
298
299

/**
 * @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.
300
 * @param verbose Print messages to stdout or not
301
302
 *
 */
303

304
void space_rebuild(struct space *s, double cell_max, int verbose) {
305

306
  int j, k, cdim[3], nr_parts = s->nr_parts, nr_gparts = s->nr_gparts;
307
  struct cell *restrict c, *restrict cells;
308
  struct part *restrict p;
309
  int *ind;
310
311
312
313
314
315
316
  double ih[3], dim[3];
  // ticks tic;

  /* 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
  cells = s->cells;

  /* Run through the particles and get their cell index. */
  // tic = getticks();
322
323
  const int ind_size = s->size_parts;
  if ((ind = (int *)malloc(sizeof(int) * ind_size)) == NULL)
324
325
326
327
328
329
330
331
332
333
    error("Failed to allocate temporary particle indices.");
  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
  for (k = 0; k < nr_parts; k++) {
335
    p = &s->parts[k];
336
337
338
339
340
    for (j = 0; j < 3; j++)
      if (p->x[j] < 0.0)
        p->x[j] += dim[j];
      else if (p->x[j] >= dim[j])
        p->x[j] -= dim[j];
341
    ind[k] =
342
        cell_getid(cdim, p->x[0] * ih[0], p->x[1] * ih[1], p->x[2] * ih[2]);
343
    cells[ind[k]].count++;
344
345
346
347
348
349
350
  }
// message( "getting particle indices took %.3f ms." , (double)(getticks() -
// tic) / CPU_TPS * 1000 );

#ifdef WITH_MPI
  /* Move non-local parts to the end of the list. */
  int nodeID = s->e->nodeID;
351
352
353
354
  for (k = 0; k < nr_parts; k++)
    if (cells[ind[k]].nodeID != nodeID) {
      cells[ind[k]].count -= 1;
      nr_parts -= 1;
355
356
357
358
359
360
      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;
361
362
363
      int t = ind[k];
      ind[k] = ind[nr_parts];
      ind[nr_parts] = t;
364
365
    }

366
367
368
  /* Exchange the strays, note that this potentially re-allocates
     the parts arrays. */
  s->nr_parts =
369
370
      nr_parts + engine_exchange_strays(s->e, nr_parts, &ind[nr_parts],
                                        s->nr_parts - nr_parts);
371
372

  /* Re-allocate the index array if needed.. */
373
374
375
  if (s->nr_parts > ind_size) {
    int *ind_new;
    if ((ind_new = (int *)malloc(sizeof(int) * s->nr_parts)) == NULL)
376
      error("Failed to allocate temporary particle indices.");
377
378
379
    memcpy(ind_new, ind, sizeof(int) * nr_parts);
    free(ind);
    ind = ind_new;
380
381
382
  }

  /* Assign each particle to its cell. */
383
  for (k = nr_parts; k < s->nr_parts; k++) {
384
    p = &s->parts[k];
385
    ind[k] =
386
        cell_getid(cdim, p->x[0] * ih[0], p->x[1] * ih[1], p->x[2] * ih[2]);
387
388
389
390
    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 ); */
391
  }
392
  nr_parts = s->nr_parts;
393
394
395
396
#endif

  /* Sort the parts according to their cells. */
  // tic = getticks();
397
  space_parts_sort(s, ind, nr_parts, 0, s->nr_cells - 1);
398
399
400
401
  // message( "parts_sort took %.3f ms." , (double)(getticks() - tic) / CPU_TPS
  // * 1000 );

  /* Re-link the gparts. */
402
  for (k = 0; k < nr_parts; k++)
403
    if (s->parts[k].gpart != NULL) s->parts[k].gpart->part = &s->parts[k];
404

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

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

  /* Run through the gravity particles and get their cell index. */
  // tic = getticks();
420
  if ((ind = (int *)malloc(sizeof(int) * s->size_gparts)) == NULL)
421
422
    error("Failed to allocate temporary particle indices.");
  for (k = 0; k < nr_gparts; k++) {
423
    struct gpart *gp = &s->gparts[k];
424
425
426
427
428
    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];
429
    ind[k] =
430
        cell_getid(cdim, gp->x[0] * ih[0], gp->x[1] * ih[1], gp->x[2] * ih[2]);
431
    cells[ind[k]].gcount++;
432
433
434
435
436
437
438
439
  }
  // 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();
440
  gparts_sort(s->gparts, ind, nr_gparts, 0, s->nr_cells - 1);
441
442
443
444
445
  // message( "gparts_sort took %.3f ms." , (double)(getticks() - tic) / CPU_TPS
  // * 1000 );

  /* Re-link the parts. */
  for (k = 0; k < nr_gparts; k++)
446
    if (s->gparts[k].id > 0) s->gparts[k].part->gpart = &s->gparts[k];
447
448

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

  /* Hook the cells up to the parts. */
  // tic = getticks();
453
454
455
  struct part *finger = s->parts;
  struct xpart *xfinger = s->xparts;
  struct gpart *gfinger = s->gparts;
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
  for (k = 0; k < s->nr_cells; k++) {
    c = &cells[k];
    c->parts = finger;
    c->xparts = xfinger;
    c->gparts = gfinger;
    finger = &finger[c->count];
    xfinger = &xfinger[c->count];
    gfinger = &gfinger[c->gcount];
  }
  // message( "hooking up cells took %.3f ms." , (double)(getticks() - tic) /
  // CPU_TPS * 1000 );

  /* At this point, we have the upper-level cells, old or new. Now make
     sure that the parts in each cell are ok. */
  // tic = getticks();
471
472
  // for (k = 0; k < s->nr_cells; k++) space_split(s, &cells[k]);
  for (k = 0; k < s->nr_cells; k++)
473
474
    scheduler_addtask(&s->e->sched, task_type_split_cell, task_subtype_none, k,
                      0, &cells[k], NULL, 0);
475
  engine_launch(s->e, s->e->nr_threads, 1 << task_type_split_cell, 0);
476

477
478
479
  // message( "space_split took %.3f ms." , (double)(getticks() - tic) / CPU_TPS
  // * 1000 );
}
480

481
/**
482
483
 * @brief Sort the particles and condensed particles according to the given
 *indices.
484
 *
Matthieu Schaller's avatar
Matthieu Schaller committed
485
 * @param s The #space.
486
487
488
489
490
 * @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.
 */
491

492
void space_parts_sort(struct space *s, int *ind, int N, int min, int max) {
493
  // Populate the global parallel_sort structure with the input data.
494
495
496
  space_sort_struct.parts = s->parts;
  space_sort_struct.xparts = s->xparts;
  space_sort_struct.ind = ind;
497
  space_sort_struct.stack_size = 2 * (max - min + 1) + 10 + s->e->nr_threads;
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
  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;

  // Add the first interval.
  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;

  // Launch the sorting tasks.
515
  engine_launch(s->e, s->e->nr_threads, (1 << task_type_psort), 0);
516
517

  /* Verify space_sort_struct. */
518
  /* for (int i = 1; i < N; i++)
519
    if (ind[i - 1] > ind[i])
520
521
      error("Sorting failed (ind[%i]=%i,ind[%i]=%i), min=%i, max=%i.", i - 1,
  ind[i - 1], i,
522
523
            ind[i], min, max);
  message("Sorting succeeded."); */
524
525
526
527

  // Clean up.
  free(space_sort_struct.stack);
}
528

529
void space_do_parts_sort() {
530

531
532
533
534
  /* Pointers to the sorting data. */
  int *ind = space_sort_struct.ind;
  struct part *parts = space_sort_struct.parts;
  struct xpart *xparts = space_sort_struct.xparts;
535

536
  /* Main loop. */
537
  while (space_sort_struct.waiting) {
538

539
    /* Grab an interval off the queue. */
540
541
    int qid =
        atomic_inc(&space_sort_struct.first) % space_sort_struct.stack_size;
542

543
    /* Wait for the entry to be ready, or for the sorting do be done. */
544
545
    while (!space_sort_struct.stack[qid].ready)
      if (!space_sort_struct.waiting) return;
546

547
    /* Get the stack entry. */
548
549
550
551
    int i = space_sort_struct.stack[qid].i;
    int j = space_sort_struct.stack[qid].j;
    int min = space_sort_struct.stack[qid].min;
    int max = space_sort_struct.stack[qid].max;
552
    space_sort_struct.stack[qid].ready = 0;
553

554
555
    /* Loop over sub-intervals. */
    while (1) {
556

557
      /* Bring beer. */
558
      const int pivot = (min + max) / 2;
559
560
      /* message("Working on interval [%i,%i] with min=%i, max=%i, pivot=%i.",
              i, j, min, max, pivot); */
561
562

      /* One pass of QuickSort's partitioning. */
563
564
      int ii = i;
      int jj = j;
565
566
567
568
      while (ii < jj) {
        while (ii <= j && ind[ii] <= pivot) ii++;
        while (jj >= i && ind[jj] > pivot) jj--;
        if (ii < jj) {
569
          int temp_i = ind[ii];
570
571
          ind[ii] = ind[jj];
          ind[jj] = temp_i;
572
          struct part temp_p = parts[ii];
573
574
          parts[ii] = parts[jj];
          parts[jj] = temp_p;
575
          struct xpart temp_xp = xparts[ii];
576
577
578
579
          xparts[ii] = xparts[jj];
          xparts[jj] = temp_xp;
        }
      }
580

581
      /* Verify space_sort_struct. */
582
583
584
585
586
587
588
589
590
591
592
593
      /* 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).");
        } */
594
595
596
597
598
599

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

        /* Recurse on the left? */
        if (jj > i && pivot > min) {
600
601
          qid = atomic_inc(&space_sort_struct.last) %
                space_sort_struct.stack_size;
602
603
          while (space_sort_struct.stack[qid].ready)
            ;
604
605
606
607
608
609
610
          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.");
611
          space_sort_struct.stack[qid].ready = 1;
612
        }
613

614
615
616
617
618
619
620
621
622
623
        /* Recurse on the right? */
        if (jj + 1 < j && pivot + 1 < max) {
          i = jj + 1;
          min = pivot + 1;
        } else
          break;

      } else {

        /* Recurse on the right? */
624
        if (pivot + 1 < max) {
625
626
          qid = atomic_inc(&space_sort_struct.last) %
                space_sort_struct.stack_size;
627
628
          while (space_sort_struct.stack[qid].ready)
            ;
629
630
631
632
633
634
635
          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.");
636
          space_sort_struct.stack[qid].ready = 1;
637
        }
638

639
640
641
642
643
644
645
        /* Recurse on the left? */
        if (jj > i && pivot > min) {
          j = jj;
          max = pivot;
        } else
          break;
      }
646

647
648
    } /* loop over sub-intervals. */

649
    atomic_dec(&space_sort_struct.waiting);
650
651

  } /* main loop. */
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
}

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;

687
688
  /* Main loop. */
  while (waiting > 0) {
689

690
691
    /* Grab an interval off the queue. */
    qid = (first++) % qstack_size;
692

693
694
695
696
697
698
699
700
701
    /* 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) {
702

703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
      /* 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;
        }
      }
721

722
      /* Verify space_sort_struct. */
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
      /* 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.");
        }
749

750
751
752
753
754
755
        /* Recurse on the right? */
        if (jj + 1 < j && pivot + 1 < max) {
          i = jj + 1;
          min = pivot + 1;
        } else
          break;
756

757
758
759
      } else {

        /* Recurse on the right? */
760
        if (pivot + 1 < max) {
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
          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. */
783

784
  /* Verify space_sort_struct. */
785
786
787
788
789
790
791
792
  /* 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);
}
793

Pedro Gonnet's avatar
Pedro Gonnet committed
794
/**
795
 * @brief Mapping function to free the sorted indices buffers.
Pedro Gonnet's avatar
Pedro Gonnet committed
796
797
 */

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

800
801
802
803
804
  if (c->sort != NULL) {
    free(c->sort);
    c->sort = NULL;
  }
}
Pedro Gonnet's avatar
Pedro Gonnet committed
805

806
807
808
/**
 * @brief Map a function to all particles in a cell recursively.
 *
809
 * @param c The #cell we are working in.
810
811
812
813
 * @param fun Function pointer to apply on the cells.
 * @param data Data passed to the function fun.
 */

Pedro Gonnet's avatar
Pedro Gonnet committed
814
815
816
817
static void rec_map_parts(struct cell *c,
                          void (*fun)(struct part *p, struct cell *c,
                                      void *data),
                          void *data) {
818
819
820
821
822
823

  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
824

825
826
827
828
829
830
  /* 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
831
/**
832
 * @brief Map a function to all particles in a space.
Pedro Gonnet's avatar
Pedro Gonnet committed
833
834
 *
 * @param s The #space we are working in.
835
836
 * @param fun Function pointer to apply on the cells.
 * @param data Data passed to the function fun.
Pedro Gonnet's avatar
Pedro Gonnet committed
837
838
 */

839
840
841
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
842

843
844
  int cid = 0;

845
  /* Call the recursive function on all higher-level cells. */
Pedro Gonnet's avatar
Pedro Gonnet committed
846
847
  for (cid = 0; cid < s->nr_cells; cid++)
    rec_map_parts(&s->cells[cid], fun, data);
848
}
849

850
851
852
853
854
855
856
857
/**
 * @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,
858
859
                                 void (*fun)(struct part *p, struct xpart *xp,
                                             struct cell *c)) {
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880

  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,
881
882
                            void (*fun)(struct part *p, struct xpart *xp,
                                        struct cell *c)) {
883
884
885
886
887
888
889
890

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

891
892
893
/**
 * @brief Map a function to all particles in a cell recursively.
 *
894
 * @param c The #cell we are working in.
895
896
897
898
 * @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.
 */
899

Pedro Gonnet's avatar
Pedro Gonnet committed
900
901
902
static void rec_map_cells_post(struct cell *c, int full,
                               void (*fun)(struct cell *c, void *data),
                               void *data) {
903

904
905
906
907
908
  int k;

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

912
913
  /* No progeny? */
  if (full || !c->split) fun(c, data);
914
}
Pedro Gonnet's avatar
Pedro Gonnet committed
915
916

/**
917
 * @brief Map a function to all particles in a aspace.
Pedro Gonnet's avatar
Pedro Gonnet committed
918
919
 *
 * @param s The #space we are working in.
920
 * @param full Map to all cells, including cells with sub-cells.
921
922
 * @param fun Function pointer to apply on the cells.
 * @param data Data passed to the function fun.
Pedro Gonnet's avatar
Pedro Gonnet committed
923
 */
924

925
void space_map_cells_post(struct space *s, int full,
Pedro Gonnet's avatar
Pedro Gonnet committed
926
                          void (*fun)(struct cell *c, void *data), void *data) {
927

928
  int cid = 0;
929

930
  /* Call the recursive function on all higher-level cells. */
Pedro Gonnet's avatar
Pedro Gonnet committed
931
932
  for (cid = 0; cid < s->nr_cells; cid++)
    rec_map_cells_post(&s->cells[cid], full, fun, data);
933
}
934

Pedro Gonnet's avatar
Pedro Gonnet committed
935
936
937
static void rec_map_cells_pre(struct cell *c, int full,
                              void (*fun)(struct cell *c, void *data),
                              void *data) {
938

939
  int k;
Pedro Gonnet's avatar
Pedro Gonnet committed
940

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

944
945
946
  /* Recurse. */
  if (c->split)
    for (k = 0; k < 8; k++)
Pedro Gonnet's avatar
Pedro Gonnet committed
947
948
      if (c->progeny[k] != NULL)
        rec_map_cells_pre(c->progeny[k], full, fun, data);
949
}
Pedro Gonnet's avatar
Pedro Gonnet committed
950

951
952
953
954
955
956
957
958
/**
 * @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
 */
959
void space_map_cells_pre(struct space *s, int full,
Pedro Gonnet's avatar
Pedro Gonnet committed
960
                         void (*fun)(struct cell *c, void *data), void *data) {
961

962
  int cid = 0;
963
964

  /* Call the recursive function on all higher-level cells. */
Pedro Gonnet's avatar
Pedro Gonnet committed
965
966
  for (cid = 0; cid < s->nr_cells; cid++)
    rec_map_cells_pre(&s->cells[cid], full, fun, data);
967
}
Pedro Gonnet's avatar
Pedro Gonnet committed
968
969
970
971
972
973
974

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

976
977
978
void space_split(struct space *s, struct cell *c) {

  int k, count = c->count, gcount = c->gcount, maxdepth = 0;
Matthieu Schaller's avatar
Matthieu Schaller committed
979
  float h, h_max = 0.0f, t_end_min = FLT_MAX, t_end_max = 0., t_end;
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
  struct cell *temp;
  struct part *p, *parts = c->parts;
  struct xpart *xp, *xparts = c->xparts;

  /* 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. */
    for (k = 0; k < 8; k++) {
      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. */
    cell_split(c);

    /* Remove any progeny with zero parts. */
    for (k = 0; k < 8; k++)
      if (c->progeny[k]->count == 0 && c->progeny[k]->gcount == 0) {
        space_recycle(s, c->progeny[k]);
        c->progeny[k] = NULL;
      } else {
        space_split(s, c->progeny[k]);
        h_max = fmaxf(h_max, c->progeny[k]->h_max);
Matthieu Schaller's avatar
Matthieu Schaller committed
1028
1029
        t_end_min = fminf(t_end_min, c->progeny[k]->t_end_min);
        t_end_max = fmaxf(t_end_max, c->progeny[k]->t_end_max);
1030
1031
1032
1033
1034
1035
        if (c->progeny[k]->maxdepth > maxdepth)
          maxdepth = c->progeny[k]->maxdepth;
      }

    /* Set the values for this cell. */
    c->h_max = h_max;
Matthieu Schaller's avatar
Matthieu Schaller committed
1036
1037
    c->t_end_min = t_end_min;
    c->t_end_max = t_end_max;
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
    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. */

    for (k = 0; k < count; k++) {
      p = &parts[k];
      xp = &xparts[k];
      xp->x_old[0] = p->x[0];
      xp->x_old[1] = p->x[1];
      xp->x_old[2] = p->x[2];
      h = p->h;
Matthieu Schaller's avatar
Matthieu Schaller committed
1059
      t_end = p->t_end;
1060
      if (h > h_max) h_max = h;
Matthieu Schaller's avatar
Matthieu Schaller committed
1061
1062
      if (t_end < t_end_min) t_end_min = t_end;
      if (t_end > t_end_max) t_end_max = t_end;
1063
    }
1064
    c->h_max = h_max;
Matthieu Schaller's avatar
Matthieu Schaller committed
1065
1066
    c->t_end_min = t_end_min;
    c->t_end_max = t_end_max;
1067
  }
1068

1069
  /* Set ownership according to the start of the parts array. */
1070
1071
  c->owner = ((c->parts - s->parts) % s->nr_parts) * s->nr_queues / s->nr_parts;
}
1072

Pedro Gonnet's avatar
Pedro Gonnet committed
1073
1074
1075
1076
1077
1078
1079
/**
 * @brief Return a used cell to the cell buffer.
 *
 * @param s The #space.
 * @param c The #cell.
 */

1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
void space_recycle(struct space *s, struct cell *c) {

  /* Lock the space. */
  lock_lock(&s->lock);

  /* Clear the cell. */
  if (lock_destroy(&c->lock) != 0) error("Failed to destroy spinlock.");

  /* Clear this cell's sort arrays. */
  if (c->sort != NULL) free(c->sort);

  /* Clear the cell data. */
  bzero(c, sizeof(struct cell));

  /* Hook this cell into the buffer. */
  c->next = s->cells_new;
  s->cells_new = c;
  s->tot_cells -= 1;

  /* Unlock the space. */
  lock_unlock_blind(&s->lock);
}
Pedro Gonnet's avatar
Pedro Gonnet committed
1102
1103
1104
1105
1106
1107
1108

/**
 * @brief Get a new empty cell.
 *
 * @param s The #space.
 */

1109
struct cell *space_getcell(struct space *s) {
Pedro Gonnet's avatar
Pedro Gonnet committed
1110

1111
1112
  struct cell *c;
  int k;
Pedro Gonnet's avatar
Pedro Gonnet committed
1113

1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
  /* Lock the space. */
  lock_lock(&s->lock);

  /* Is the buffer empty? */
  if (s->cells_new == NULL) {
    if (posix_memalign((void *)&s->cells_new, 64,
                       space_cellallocchunk * sizeof(struct cell)) != 0)
      error("Failed to allocate more cells.");
    bzero(s->cells_new, space_cellallocchunk * sizeof(struct cell));
    for (k = 0; k < space_cellallocchunk - 1; k++)
      s->cells_new[k].next = &s->cells_new[k + 1];
    s->cells_new[space_cellallocchunk - 1].next = NULL;
  }

  /* Pick off the next cell. */
  c = s->cells_new;
  s->cells_new = c->next;
  s->tot_cells += 1;
Pedro Gonnet's avatar
Pedro Gonnet committed
1132

1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
  /* Unlock the space. */
  lock_unlock_blind(&s->lock);

  /* Init some things in the cell. */
  bzero(c, sizeof(struct cell));
  c->nodeID = -1;
  if (lock_init(&c->lock) != 0 || lock_init(&c->glock) != 0)
    error("Failed to initialize cell spinlocks.");

  return c;
}
Pedro Gonnet's avatar
Pedro Gonnet committed
1144
1145
1146
1147

/**
 * @brief Split the space into cells given the array of particles.
 *
1148
 * @param s The #space to initialize.
Pedro Gonnet's avatar
Pedro Gonnet committed
1149
1150
1151
1152
 * @param dim Spatial dimensions of the domain.
 * @param parts Pointer to an array of #part.
 * @param N The number of parts in the space.
 * @param periodic flag whether the domain is periodic or not.
1153
 * @param h_max The maximal interaction radius.
1154
 * @param verbose Print messages to stdout or not
Pedro Gonnet's avatar
Pedro Gonnet committed
1155
1156
 *
 * Makes a grid of edge length > r_max and fills the particles
1157
 * into the respective cells. Cells containing more than #space_splitsize
Pedro Gonnet's avatar
Pedro Gonnet committed
1158
1159
1160
1161
 * parts with a cutoff below half the cell width are then split
 * recursively.
 */

1162
void space_init(struct space *s, double dim[3], struct part *parts, int N,
1163
                int periodic, double h_max, int verbose) {
1164

1165
  /* Store everything in the space. */
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
  s->dim[0] = dim[0];
  s->dim[1] = dim[1];
  s->dim[2] = dim[2];
  s->periodic = periodic;
  s->nr_parts = N;
  s->size_parts = N;
  s->parts = parts;
  s->cell_min = h_max;
  s->nr_queues = 1;
  s->size_parts_foreign = 0;

  /* Check that all the particle positions are reasonable, wrap if periodic. */
  if (periodic) {
    for (int k = 0; k < N; k++)
      for (int j = 0; j < 3; j++) {
        while (parts[k].x[j] < 0) parts[k].x[j] += dim[j];
        while (parts[k].x[j] >= dim[j]) parts[k].x[j] -= dim[j];
1183
      }
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
  } else {
    for (int k = 0; k < N; k++)
      for (int j = 0; j < 3; j++)
        if (parts[k].x[j] < 0 || parts[k].x[j] >= dim[j])
          error("Not all particles are within the specified domain.");
  }

  /* Allocate the xtra parts array. */
  if (posix_memalign((void *)&s->xparts, part_align,
                     N * sizeof(struct xpart)) != 0)
    error("Failed to allocate xparts.");
  bzero(s->xparts, N * sizeof(struct xpart));

  /* For now, clone the parts to make gparts. */
  if (posix_memalign((void *)&s->gparts, part_align,
                     N * sizeof(struct gpart)) != 0)
    error("Failed to allocate gparts.");
  bzero(s->gparts, N * sizeof(struct gpart));
  /* for ( int k = 0 ; k < N ; k++ ) {
      s->gparts[k].x[0] = s->parts[k].x[0];
      s->gparts[k].x[1] = s->parts[k].x[1];
      s->gparts[k].x[2] = s->parts[k].x[2];
      s->gparts[k].v[0] = s->parts[k].v[0];
      s->gparts[k].v[1] = s->parts[k].v[1];
      s->gparts[k].v[2] = s->parts[k].v[2];
      s->gparts[k].mass = s->parts[k].mass;
      s->gparts[k].dt = s->parts[k].dt;
      s->gparts[k].id = s->parts[k].id;
      s->gparts[k].part = &s->parts[k];
      s->parts[k].gpart = &s->gparts[k];
1214
      }
1215
1216
1217
1218
1219
1220
  s->nr_gparts = s->nr_parts; */
  s->nr_gparts = 0;
  s->size_gparts = s->size_parts;

  /* Init the space lock. */
  if (lock_init(&s->lock) != 0) error("Failed to create space spin-lock.");
Pedro Gonnet's avatar
Pedro Gonnet committed
1221

1222
  /* Build the cells and the tasks. */
1223
  space_regrid(s, h_max, verbose);
1224
}
1225
1226
1227
1228
1229
1230
1231
1232

/**
 * @brief Cleans-up all the cell links in the space
 *
 * Expensive funtion. Should only be used for debugging purposes.
 */
void space_link_cleanup(struct space *s) {

Matthieu Schaller's avatar
Matthieu Schaller committed
1233
  /* Recursively apply the cell link cleaning routine */
1234
1235
  space_map_cells_pre(s, 1, cell_clean_links, NULL);
}