task.c 39.2 KB
Newer Older
1
/*******************************************************************************
2
 * This file is part of SWIFT.
3
 * Copyright (c) 2012 Pedro Gonnet (pedro.gonnet@durham.ac.uk)
4
5
6
7
 *                    Matthieu Schaller (matthieu.schaller@durham.ac.uk)
 *               2015 Peter W. Draper (p.w.draper@durham.ac.uk)
 *               2016 John A. Regan (john.a.regan@durham.ac.uk)
 *                    Tom Theuns (tom.theuns@durham.ac.uk)
8
 *
9
10
11
12
 * 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.
13
 *
14
15
16
17
 * 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.
18
 *
19
20
 * 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/>.
21
 *
22
23
24
25
26
27
28
29
30
 ******************************************************************************/

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

/* Some standard headers. */
#include <float.h>
#include <limits.h>
#include <sched.h>
31
32
33
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
34

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

40
41
42
/* This object's header. */
#include "task.h"

43
/* Local headers. */
Pedro Gonnet's avatar
Pedro Gonnet committed
44
#include "atomic.h"
45
#include "engine.h"
46
#include "error.h"
47
#include "inline.h"
48
#include "lock.h"
Peter W. Draper's avatar
Peter W. Draper committed
49
#include "mpiuse.h"
50
51

/* Task type names. */
52
53
54
55
56
57
58
59
60
61
62
63
64
const char *taskID_names[task_type_count] = {"none",
                                             "sort",
                                             "self",
                                             "pair",
                                             "sub_self",
                                             "sub_pair",
                                             "init_grav",
                                             "init_grav_out",
                                             "ghost_in",
                                             "ghost",
                                             "ghost_out",
                                             "extra_ghost",
                                             "drift_part",
65
                                             "drift_spart",
66
                                             "drift_bpart",
67
68
                                             "drift_gpart",
                                             "drift_gpart_out",
69
                                             "end_hydro_force",
70
71
72
                                             "kick1",
                                             "kick2",
                                             "timestep",
73
                                             "timestep_limiter",
74
					     "timestep_sync",
75
76
                                             "limiter_in",
                                             "limiter_out",
77
78
79
80
81
82
83
                                             "send",
                                             "recv",
                                             "grav_long_range",
                                             "grav_mm",
                                             "grav_down_in",
                                             "grav_down",
                                             "grav_mesh",
84
                                             "grav_end_force",
85
86
                                             "cooling",
                                             "star_formation",
87
88
                                             "star_formation_in",
                                             "star_formation_out",
89
                                             "logger",
90
91
                                             "stars_in",
                                             "stars_out",
92
93
                                             "stars_ghost_in",
                                             "stars_ghost",
James Willis's avatar
James Willis committed
94
                                             "stars_ghost_out",
95
                                             "stars_sort",
96
                                             "stars_resort",
97
98
                                             "bh_in",
                                             "bh_out",
99
100
101
                                             "bh_density_ghost",
                                             "bh_swallow_ghost1",
                                             "bh_swallow_ghost2",
102
                                             "bh_swallow_ghost3",
103
                                             "fof_self",
James Willis's avatar
James Willis committed
104
                                             "fof_pair"};
105

106
/* Sub-task type names. */
107
const char *subtaskID_names[task_subtype_count] = {
108
109
110
111
112
113
114
115
    "none",       "density",      "gradient",       "force",
    "limiter",    "grav",         "external_grav",  "tend_part",
    "tend_gpart", "tend_spart",   "tend_bpart",     "xv",
    "rho",        "part_swallow", "bpart_merger",   "gpart",
    "multipole",  "spart",        "stars_density",  "stars_feedback",
    "sf_count",   "bpart_rho",    "bpart_swallow",  "bpart_feedback",
    "bh_density", "bh_swallow",   "do_gas_swallow", "do_bh_swallow",
    "bh_feedback"};
116

117
118
119
120
121
#ifdef WITH_MPI
/* MPI communicators for the subtypes. */
MPI_Comm subtaskMPI_comms[task_subtype_count];
#endif

122
123
/**
 * @brief Computes the overlap between the parts array of two given cells.
124
 *
Matthieu Schaller's avatar
Matthieu Schaller committed
125
126
127
 * @param TYPE is the type of parts (e.g. #part, #gpart, #spart)
 * @param ARRAY is the array of this specific type.
 * @param COUNT is the number of elements in the array.
128
 */
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
#define TASK_CELL_OVERLAP(TYPE, ARRAY, COUNT)                               \
  __attribute__((always_inline))                                            \
      INLINE static size_t task_cell_overlap_##TYPE(                        \
          const struct cell *restrict ci, const struct cell *restrict cj) { \
                                                                            \
    if (ci == NULL || cj == NULL) return 0;                                 \
                                                                            \
    if (ci->ARRAY <= cj->ARRAY &&                                           \
        ci->ARRAY + ci->COUNT >= cj->ARRAY + cj->COUNT) {                   \
      return cj->COUNT;                                                     \
    } else if (cj->ARRAY <= ci->ARRAY &&                                    \
               cj->ARRAY + cj->COUNT >= ci->ARRAY + ci->COUNT) {            \
      return ci->COUNT;                                                     \
    }                                                                       \
                                                                            \
    return 0;                                                               \
  }
146

147
TASK_CELL_OVERLAP(part, hydro.parts, hydro.count);
148
149
TASK_CELL_OVERLAP(gpart, grav.parts, grav.count);
TASK_CELL_OVERLAP(spart, stars.parts, stars.count);
150
TASK_CELL_OVERLAP(bpart, black_holes.parts, black_holes.count);
Loic Hausammann's avatar
Loic Hausammann committed
151

152
153
154
155
156
/**
 * @brief Returns the #task_actions for a given task.
 *
 * @param t The #task.
 */
157
158
__attribute__((always_inline)) INLINE static enum task_actions task_acts_on(
    const struct task *t) {
159
160
161
162
163
164
165

  switch (t->type) {

    case task_type_none:
      return task_action_none;
      break;

166
    case task_type_drift_part:
167
168
    case task_type_sort:
    case task_type_ghost:
169
    case task_type_extra_ghost:
170
    case task_type_timestep_limiter:
Stefan Arridge's avatar
Stefan Arridge committed
171
    case task_type_cooling:
172
    case task_type_end_hydro_force:
173
174
175
      return task_action_part;
      break;

176
177
178
    case task_type_star_formation:
      return task_action_all;

179
    case task_type_drift_spart:
180
    case task_type_stars_ghost:
Loic Hausammann's avatar
Loic Hausammann committed
181
    case task_type_stars_sort:
182
    case task_type_stars_resort:
Loic Hausammann's avatar
Loic Hausammann committed
183
184
185
      return task_action_spart;
      break;

186
    case task_type_drift_bpart:
187
    case task_type_bh_density_ghost:
188
    case task_type_bh_swallow_ghost3:
189
190
191
      return task_action_bpart;
      break;

192
193
194
195
196
197
198
    case task_type_self:
    case task_type_pair:
    case task_type_sub_self:
    case task_type_sub_pair:
      switch (t->subtype) {

        case task_subtype_density:
199
        case task_subtype_gradient:
200
        case task_subtype_force:
201
        case task_subtype_limiter:
202
203
204
          return task_action_part;
          break;

205
        case task_subtype_stars_density:
Alexei Borissov's avatar
Alexei Borissov committed
206
        case task_subtype_stars_feedback:
207
208
          return task_action_all;
          break;
209

210
211
        case task_subtype_bh_density:
        case task_subtype_bh_feedback:
212
        case task_subtype_bh_swallow:
213
        case task_subtype_do_gas_swallow:
214
215
216
          return task_action_all;
          break;

217
218
219
220
        case task_subtype_do_bh_swallow:
          return task_action_bpart;
          break;

221
        case task_subtype_grav:
222
        case task_subtype_external_grav:
223
224
225
226
          return task_action_gpart;
          break;

        default:
227
228
229
230
#ifdef SWIFT_DEBUG_CHECKS
          error("Unknown task_action for task %s/%s", taskID_names[t->type],
                subtaskID_names[t->subtype]);
#endif
231
232
233
234
235
          return task_action_none;
          break;
      }
      break;

236
237
    case task_type_kick1:
    case task_type_kick2:
Loikki's avatar
Loikki committed
238
    case task_type_logger:
James Willis's avatar
James Willis committed
239
240
    case task_type_fof_self:
    case task_type_fof_pair:
241
    case task_type_timestep:
242
243
    case task_type_send:
    case task_type_recv:
244
      if (t->ci->hydro.count > 0 && t->ci->grav.count > 0)
245
        return task_action_all;
246
      else if (t->ci->hydro.count > 0)
247
        return task_action_part;
248
      else if (t->ci->grav.count > 0)
249
        return task_action_gpart;
250
251
      else {
#ifdef SWIFT_DEBUG_CHECKS
252
        error("Task without particles");
253
254
#endif
      }
255
256
      break;

257
    case task_type_init_grav:
258
    case task_type_grav_mm:
259
    case task_type_grav_long_range:
260
261
262
      return task_action_multipole;
      break;

263
    case task_type_drift_gpart:
264
    case task_type_grav_down:
265
    case task_type_end_grav_force:
266
    case task_type_grav_mesh:
267
      return task_action_gpart;
268
      break;
269

270
    default:
271
272
273
274
#ifdef SWIFT_DEBUG_CHECKS
      error("Unknown task_action for task %s/%s", taskID_names[t->type],
            subtaskID_names[t->subtype]);
#endif
275
276
277
      return task_action_none;
      break;
  }
278

279
280
281
282
283
#ifdef SWIFT_DEBUG_CHECKS
  error("Unknown task_action for task %s/%s", taskID_names[t->type],
        subtaskID_names[t->subtype]);
#endif
  /* Silence compiler warnings. We should never get here. */
284
  return task_action_none;
285
286
}

287
288
289
290
291
292
293
/**
 * @brief Compute the Jaccard similarity of the data used by two
 *        different tasks.
 *
 * @param ta The first #task.
 * @param tb The second #task.
 */
294
295
float task_overlap(const struct task *restrict ta,
                   const struct task *restrict tb) {
296
297
298
299
300
301

  if (ta == NULL || tb == NULL) return 0.f;

  const enum task_actions ta_act = task_acts_on(ta);
  const enum task_actions tb_act = task_acts_on(tb);

302
303
  /* First check if any of the two tasks are of a type that don't
     use cells. */
304
305
306
307
308
  if (ta_act == task_action_none || tb_act == task_action_none) return 0.f;

  const int ta_part = (ta_act == task_action_part || ta_act == task_action_all);
  const int ta_gpart =
      (ta_act == task_action_gpart || ta_act == task_action_all);
309
310
  const int ta_spart =
      (ta_act == task_action_spart || ta_act == task_action_all);
311
312
  const int ta_bpart =
      (ta_act == task_action_bpart || ta_act == task_action_all);
313
314
315
  const int tb_part = (tb_act == task_action_part || tb_act == task_action_all);
  const int tb_gpart =
      (tb_act == task_action_gpart || tb_act == task_action_all);
316
317
  const int tb_spart =
      (tb_act == task_action_spart || tb_act == task_action_all);
318
319
  const int tb_bpart =
      (tb_act == task_action_bpart || tb_act == task_action_all);
320
321
322
323
324
325

  /* In the case where both tasks act on parts */
  if (ta_part && tb_part) {

    /* Compute the union of the cell data. */
    size_t size_union = 0;
326
327
328
329
    if (ta->ci != NULL) size_union += ta->ci->hydro.count;
    if (ta->cj != NULL) size_union += ta->cj->hydro.count;
    if (tb->ci != NULL) size_union += tb->ci->hydro.count;
    if (tb->cj != NULL) size_union += tb->cj->hydro.count;
330

331
    if (size_union == 0) return 0.f;
332

333
334
335
336
337
338
339
340
341
342
343
344
345
346
    /* Compute the intersection of the cell data. */
    const size_t size_intersect = task_cell_overlap_part(ta->ci, tb->ci) +
                                  task_cell_overlap_part(ta->ci, tb->cj) +
                                  task_cell_overlap_part(ta->cj, tb->ci) +
                                  task_cell_overlap_part(ta->cj, tb->cj);

    return ((float)size_intersect) / (size_union - size_intersect);
  }

  /* In the case where both tasks act on gparts */
  else if (ta_gpart && tb_gpart) {

    /* Compute the union of the cell data. */
    size_t size_union = 0;
347
348
349
350
    if (ta->ci != NULL) size_union += ta->ci->grav.count;
    if (ta->cj != NULL) size_union += ta->cj->grav.count;
    if (tb->ci != NULL) size_union += tb->ci->grav.count;
    if (tb->cj != NULL) size_union += tb->cj->grav.count;
351

352
353
    if (size_union == 0) return 0.f;

354
355
356
357
358
359
360
361
    /* Compute the intersection of the cell data. */
    const size_t size_intersect = task_cell_overlap_gpart(ta->ci, tb->ci) +
                                  task_cell_overlap_gpart(ta->ci, tb->cj) +
                                  task_cell_overlap_gpart(ta->cj, tb->ci) +
                                  task_cell_overlap_gpart(ta->cj, tb->cj);

    return ((float)size_intersect) / (size_union - size_intersect);
  }
362

Loic Hausammann's avatar
Loic Hausammann committed
363
364
365
366
367
  /* In the case where both tasks act on sparts */
  else if (ta_spart && tb_spart) {

    /* Compute the union of the cell data. */
    size_t size_union = 0;
368
369
370
371
    if (ta->ci != NULL) size_union += ta->ci->stars.count;
    if (ta->cj != NULL) size_union += ta->cj->stars.count;
    if (tb->ci != NULL) size_union += tb->ci->stars.count;
    if (tb->cj != NULL) size_union += tb->cj->stars.count;
Loic Hausammann's avatar
Loic Hausammann committed
372

373
    if (size_union == 0) return 0.f;
Loic Hausammann's avatar
Loic Hausammann committed
374

Loic Hausammann's avatar
Loic Hausammann committed
375
376
377
378
379
380
381
382
383
    /* Compute the intersection of the cell data. */
    const size_t size_intersect = task_cell_overlap_spart(ta->ci, tb->ci) +
                                  task_cell_overlap_spart(ta->ci, tb->cj) +
                                  task_cell_overlap_spart(ta->cj, tb->ci) +
                                  task_cell_overlap_spart(ta->cj, tb->cj);

    return ((float)size_intersect) / (size_union - size_intersect);
  }

384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
  /* In the case where both tasks act on bparts */
  else if (ta_bpart && tb_bpart) {

    /* Compute the union of the cell data. */
    size_t size_union = 0;
    if (ta->ci != NULL) size_union += ta->ci->black_holes.count;
    if (ta->cj != NULL) size_union += ta->cj->black_holes.count;
    if (tb->ci != NULL) size_union += tb->ci->black_holes.count;
    if (tb->cj != NULL) size_union += tb->cj->black_holes.count;

    if (size_union == 0) return 0.f;

    /* Compute the intersection of the cell data. */
    const size_t size_intersect = task_cell_overlap_bpart(ta->ci, tb->ci) +
                                  task_cell_overlap_bpart(ta->ci, tb->cj) +
                                  task_cell_overlap_bpart(ta->cj, tb->ci) +
                                  task_cell_overlap_bpart(ta->cj, tb->cj);

    return ((float)size_intersect) / (size_union - size_intersect);
  }

405
406
  /* Else, no overlap */
  return 0.f;
407
}
408

409
410
/**
 * @brief Unlock the cell held by this task.
411
 *
412
413
 * @param t The #task.
 */
414
415
void task_unlock(struct task *t) {

416
417
  const enum task_types type = t->type;
  const enum task_subtypes subtype = t->subtype;
418
419
  struct cell *ci = t->ci, *cj = t->cj;

420
  /* Act based on task type. */
421
422
  switch (type) {

423
424
    case task_type_kick1:
    case task_type_kick2:
425
    case task_type_logger:
426
    case task_type_timestep:
427
428
429
      cell_unlocktree(ci);
      cell_gunlocktree(ci);
      break;
Matthieu Schaller's avatar
Matthieu Schaller committed
430

431
    case task_type_drift_part:
432
    case task_type_sort:
433
    case task_type_ghost:
434
    case task_type_extra_ghost:
435
    case task_type_end_hydro_force:
436
    case task_type_timestep_limiter:
437
    case task_type_timestep_sync:
438
439
440
      cell_unlocktree(ci);
      break;

441
    case task_type_drift_gpart:
442
    case task_type_grav_mesh:
443
    case task_type_end_grav_force:
444
445
446
      cell_gunlocktree(ci);
      break;

Loic Hausammann's avatar
Loic Hausammann committed
447
    case task_type_stars_sort:
448
    case task_type_stars_resort:
Loic Hausammann's avatar
Loic Hausammann committed
449
450
451
      cell_sunlocktree(ci);
      break;

452
    case task_type_self:
453
    case task_type_sub_self:
454
455
      if (subtype == task_subtype_grav) {
        cell_gunlocktree(ci);
456
        cell_munlocktree(ci);
457
458
      } else if ((subtype == task_subtype_stars_density) ||
                 (subtype == task_subtype_stars_feedback)) {
Alexei Borissov's avatar
Alexei Borissov committed
459
460
        cell_sunlocktree(ci);
        cell_unlocktree(ci);
461
462
463
      } else if ((subtype == task_subtype_bh_density) ||
                 (subtype == task_subtype_bh_feedback) ||
                 (subtype == task_subtype_bh_swallow) ||
464
                 (subtype == task_subtype_do_gas_swallow)) {
465
466
        cell_bunlocktree(ci);
        cell_unlocktree(ci);
467
468
      } else if (subtype == task_subtype_do_bh_swallow) {
        cell_bunlocktree(ci);
469
470
471
      } else {
        cell_unlocktree(ci);
      }
472
      break;
473

474
    case task_type_pair:
475
    case task_type_sub_pair:
476
477
478
      if (subtype == task_subtype_grav) {
        cell_gunlocktree(ci);
        cell_gunlocktree(cj);
479
480
        cell_munlocktree(ci);
        cell_munlocktree(cj);
481
482
      } else if ((subtype == task_subtype_stars_density) ||
                 (subtype == task_subtype_stars_feedback)) {
Alexei Borissov's avatar
Alexei Borissov committed
483
484
485
486
        cell_sunlocktree(ci);
        cell_sunlocktree(cj);
        cell_unlocktree(ci);
        cell_unlocktree(cj);
487
488
489
      } else if ((subtype == task_subtype_bh_density) ||
                 (subtype == task_subtype_bh_feedback) ||
                 (subtype == task_subtype_bh_swallow) ||
490
                 (subtype == task_subtype_do_gas_swallow)) {
491
492
493
494
        cell_bunlocktree(ci);
        cell_bunlocktree(cj);
        cell_unlocktree(ci);
        cell_unlocktree(cj);
495
496
497
      } else if (subtype == task_subtype_do_bh_swallow) {
        cell_bunlocktree(ci);
        cell_bunlocktree(cj);
498
499
500
501
502
503
      } else {
        cell_unlocktree(ci);
        cell_unlocktree(cj);
      }
      break;

504
    case task_type_grav_down:
505
      cell_gunlocktree(ci);
506
507
508
      cell_munlocktree(ci);
      break;

509
    case task_type_grav_long_range:
510
      cell_munlocktree(ci);
511
      break;
512

513
514
515
516
517
    case task_type_grav_mm:
      cell_munlocktree(ci);
      cell_munlocktree(cj);
      break;

518
519
520
521
    case task_type_star_formation:
      cell_unlocktree(ci);
      cell_sunlocktree(ci);
      cell_gunlocktree(ci);
522
      break;
523

524
525
526
527
    default:
      break;
  }
}
528
529
530
531
532
533

/**
 * @brief Try to lock the cells associated with this task.
 *
 * @param t the #task.
 */
534
535
int task_lock(struct task *t) {

536
537
  const enum task_types type = t->type;
  const enum task_subtypes subtype = t->subtype;
538
  struct cell *ci = t->ci, *cj = t->cj;
539
540
541
542
#ifdef WITH_MPI
  int res = 0, err = 0;
  MPI_Status stat;
#endif
543

544
  switch (type) {
545

546
547
548
    /* Communication task? */
    case task_type_recv:
    case task_type_send:
549
#ifdef WITH_MPI
550
551
552
553
554
      /* Check the status of the MPI request. */
      if ((err = MPI_Test(&t->req, &res, &stat)) != MPI_SUCCESS) {
        char buff[MPI_MAX_ERROR_STRING];
        int len;
        MPI_Error_string(err, buff, &len);
555
556
557
558
        error(
            "Failed to test request on send/recv task (type=%s/%s tag=%lld, "
            "%s).",
            taskID_names[t->type], subtaskID_names[t->subtype], t->flags, buff);
559
      }
Peter W. Draper's avatar
Peter W. Draper committed
560
561

      /* And log deactivation, if logging enabled. */
562
563
564
      if (res) {
        mpiuse_log_allocation(t->type, t->subtype, &t->req, 0, 0, 0, 0);
      }
Peter W. Draper's avatar
Peter W. Draper committed
565

566
      return res;
567
#else
568
      error("SWIFT was not compiled with MPI support.");
569
#endif
570
      break;
571

572
573
    case task_type_kick1:
    case task_type_kick2:
Loikki's avatar
Loikki committed
574
    case task_type_logger:
575
    case task_type_timestep:
576
      if (ci->hydro.hold || ci->grav.phold) return 0;
577
578
      if (cell_locktree(ci) != 0) return 0;
      if (cell_glocktree(ci) != 0) {
Matthieu Schaller's avatar
Matthieu Schaller committed
579
580
        cell_unlocktree(ci);
        return 0;
581
582
583
      }
      break;

584
    case task_type_drift_part:
585
    case task_type_sort:
586
    case task_type_ghost:
587
    case task_type_extra_ghost:
588
    case task_type_end_hydro_force:
589
    case task_type_timestep_limiter:
590
    case task_type_timestep_sync:
591
      if (ci->hydro.hold) return 0;
592
593
      if (cell_locktree(ci) != 0) return 0;
      break;
594

Loic Hausammann's avatar
Loic Hausammann committed
595
    case task_type_stars_sort:
596
    case task_type_stars_resort:
Loic Hausammann's avatar
Loic Hausammann committed
597
598
599
600
      if (ci->stars.hold) return 0;
      if (cell_slocktree(ci) != 0) return 0;
      break;

601
    case task_type_drift_gpart:
602
    case task_type_end_grav_force:
603
    case task_type_grav_mesh:
604
      if (ci->grav.phold) return 0;
605
606
607
      if (cell_glocktree(ci) != 0) return 0;
      break;

608
    case task_type_self:
609
    case task_type_sub_self:
610
      if (subtype == task_subtype_grav) {
611
        /* Lock the gparts and the m-pole */
612
        if (ci->grav.phold || ci->grav.mhold) return 0;
613
614
615
616
617
618
        if (cell_glocktree(ci) != 0)
          return 0;
        else if (cell_mlocktree(ci) != 0) {
          cell_gunlocktree(ci);
          return 0;
        }
619
620
      } else if ((subtype == task_subtype_stars_density) ||
                 (subtype == task_subtype_stars_feedback)) {
Alexei Borissov's avatar
Alexei Borissov committed
621
622
623
624
625
626
627
        if (ci->stars.hold) return 0;
        if (ci->hydro.hold) return 0;
        if (cell_slocktree(ci) != 0) return 0;
        if (cell_locktree(ci) != 0) {
          cell_sunlocktree(ci);
          return 0;
        }
628
629
630
      } else if ((subtype == task_subtype_bh_density) ||
                 (subtype == task_subtype_bh_feedback) ||
                 (subtype == task_subtype_bh_swallow) ||
631
                 (subtype == task_subtype_do_gas_swallow)) {
632
633
634
635
636
637
638
        if (ci->black_holes.hold) return 0;
        if (ci->hydro.hold) return 0;
        if (cell_blocktree(ci) != 0) return 0;
        if (cell_locktree(ci) != 0) {
          cell_bunlocktree(ci);
          return 0;
        }
639
640
641
      } else if (subtype == task_subtype_do_bh_swallow) {
        if (ci->black_holes.hold) return 0;
        if (cell_blocktree(ci) != 0) return 0;
Alexei Borissov's avatar
Alexei Borissov committed
642
      } else { /* subtype == hydro */
Loic Hausammann's avatar
Loic Hausammann committed
643
        if (ci->hydro.hold) return 0;
644
645
646
        if (cell_locktree(ci) != 0) return 0;
      }
      break;
647

648
    case task_type_pair:
649
    case task_type_sub_pair:
650
      if (subtype == task_subtype_grav) {
651
        /* Lock the gparts and the m-pole in both cells */
652
        if (ci->grav.phold || cj->grav.phold) return 0;
653
654
655
656
        if (cell_glocktree(ci) != 0) return 0;
        if (cell_glocktree(cj) != 0) {
          cell_gunlocktree(ci);
          return 0;
657
658
659
660
661
662
663
664
665
        } else if (cell_mlocktree(ci) != 0) {
          cell_gunlocktree(ci);
          cell_gunlocktree(cj);
          return 0;
        } else if (cell_mlocktree(cj) != 0) {
          cell_gunlocktree(ci);
          cell_gunlocktree(cj);
          cell_munlocktree(ci);
          return 0;
666
        }
667
668
      } else if ((subtype == task_subtype_stars_density) ||
                 (subtype == task_subtype_stars_feedback)) {
Alexei Borissov's avatar
Alexei Borissov committed
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
        /* Lock the stars and the gas particles in both cells */
        if (ci->stars.hold || cj->stars.hold) return 0;
        if (ci->hydro.hold || cj->hydro.hold) return 0;
        if (cell_slocktree(ci) != 0) return 0;
        if (cell_slocktree(cj) != 0) {
          cell_sunlocktree(ci);
          return 0;
        }
        if (cell_locktree(ci) != 0) {
          cell_sunlocktree(ci);
          cell_sunlocktree(cj);
          return 0;
        }
        if (cell_locktree(cj) != 0) {
          cell_sunlocktree(ci);
          cell_sunlocktree(cj);
          cell_unlocktree(ci);
          return 0;
        }
688
689
690
      } else if ((subtype == task_subtype_bh_density) ||
                 (subtype == task_subtype_bh_feedback) ||
                 (subtype == task_subtype_bh_swallow) ||
691
                 (subtype == task_subtype_do_gas_swallow)) {
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
        /* Lock the BHs and the gas particles in both cells */
        if (ci->black_holes.hold || cj->black_holes.hold) return 0;
        if (ci->hydro.hold || cj->hydro.hold) return 0;
        if (cell_blocktree(ci) != 0) return 0;
        if (cell_blocktree(cj) != 0) {
          cell_bunlocktree(ci);
          return 0;
        }
        if (cell_locktree(ci) != 0) {
          cell_bunlocktree(ci);
          cell_bunlocktree(cj);
          return 0;
        }
        if (cell_locktree(cj) != 0) {
          cell_bunlocktree(ci);
          cell_bunlocktree(cj);
          cell_unlocktree(ci);
          return 0;
        }
711
712
713
714
715
716
717
      } else if (subtype == task_subtype_do_bh_swallow) {
        if (ci->black_holes.hold || cj->black_holes.hold) return 0;
        if (cell_blocktree(ci) != 0) return 0;
        if (cell_blocktree(cj) != 0) {
          cell_bunlocktree(ci);
          return 0;
        }
Alexei Borissov's avatar
Alexei Borissov committed
718
      } else { /* subtype == hydro */
719
        /* Lock the parts in both cells */
720
        if (ci->hydro.hold || cj->hydro.hold) return 0;
721
722
723
724
725
726
727
        if (cell_locktree(ci) != 0) return 0;
        if (cell_locktree(cj) != 0) {
          cell_unlocktree(ci);
          return 0;
        }
      }
      break;
728

729
730
    case task_type_grav_down:
      /* Lock the gparts and the m-poles */
731
      if (ci->grav.phold || ci->grav.mhold) return 0;
732
733
734
735
736
737
738
739
      if (cell_glocktree(ci) != 0)
        return 0;
      else if (cell_mlocktree(ci) != 0) {
        cell_gunlocktree(ci);
        return 0;
      }
      break;

740
    case task_type_grav_long_range:
741
      /* Lock the m-poles */
742
      if (ci->grav.mhold) return 0;
743
      if (cell_mlocktree(ci) != 0) return 0;
Matthieu Schaller's avatar
Matthieu Schaller committed
744
745
      break;

746
747
    case task_type_grav_mm:
      /* Lock both m-poles */
748
      if (ci->grav.mhold || cj->grav.mhold) return 0;
749
750
751
752
753
      if (cell_mlocktree(ci) != 0) return 0;
      if (cell_mlocktree(cj) != 0) {
        cell_munlocktree(ci);
        return 0;
      }
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
      break;

    case task_type_star_formation:
      /* Lock the gas, gravity and star particles */
      if (ci->hydro.hold || ci->stars.hold || ci->grav.phold) return 0;
      if (cell_locktree(ci) != 0) return 0;
      if (cell_slocktree(ci) != 0) {
        cell_unlocktree(ci);
        return 0;
      }
      if (cell_glocktree(ci) != 0) {
        cell_unlocktree(ci);
        cell_sunlocktree(ci);
        return 0;
      }
769

770
771
    default:
      break;
772
773
774
775
776
  }

  /* If we made it this far, we've got a lock. */
  return 1;
}
777

778
779
780
781
782
783
784
785
786
787
788
/**
 * @brief Print basic information about a task.
 *
 * @param t The #task.
 */
void task_print(const struct task *t) {

  message("Type:'%s' sub_type:'%s' wait=%d nr_unlocks=%d skip=%d",
          taskID_names[t->type], subtaskID_names[t->subtype], t->wait,
          t->nr_unlock_tasks, t->skip);
}
789

790
791
792
793
794
795
/**
 * @brief Get the group name of a task.
 *
 * This is used to group tasks with similar actions in the task dependency
 * graph.
 *
796
 * @param type The #task type.
797
 * @param subtype The #task subtype.
798
 * @param cluster (return) The group name (should be allocated)
799
 */
800
void task_get_group_name(int type, int subtype, char *cluster) {
801

802
803
  if (type == task_type_grav_long_range || type == task_type_grav_mm ||
      type == task_type_grav_mesh) {
804
805
806
807
808

    strcpy(cluster, "Gravity");
    return;
  }

809
  switch (subtype) {
810
811
812
813
    case task_subtype_density:
      strcpy(cluster, "Density");
      break;
    case task_subtype_gradient:
814
      if (type == task_type_send || type == task_type_recv) {
815
816
817
        strcpy(cluster, "None");
      } else {
        strcpy(cluster, "Gradient");
818
      }
819
820
821
822
823
824
825
      break;
    case task_subtype_force:
      strcpy(cluster, "Force");
      break;
    case task_subtype_grav:
      strcpy(cluster, "Gravity");
      break;
826
827
828
    case task_subtype_limiter:
      strcpy(cluster, "Timestep_limiter");
      break;
829
    case task_subtype_stars_density:
830
831
832
833
      strcpy(cluster, "StarsDensity");
      break;
    case task_subtype_stars_feedback:
      strcpy(cluster, "StarsFeedback");
834
      break;
835
836
837
    case task_subtype_bh_density:
      strcpy(cluster, "BHDensity");
      break;
838
839
840
    case task_subtype_bh_swallow:
      strcpy(cluster, "BHSwallow");
      break;
841
842
843
844
845
    case task_subtype_do_gas_swallow:
      strcpy(cluster, "DoGasSwallow");
      break;
    case task_subtype_do_bh_swallow:
      strcpy(cluster, "DoBHSwallow");
846
      break;
847
848
849
    case task_subtype_bh_feedback:
      strcpy(cluster, "BHFeedback");
      break;
850
851
852
853
854
855
856
857
858
859
860
861
862
    default:
      strcpy(cluster, "None");
      break;
  }
}

/**
 * @brief Generate the full name of a #task.
 *
 * @param type The #task type.
 * @param subtype The #task type.
 * @param name (return) The formatted string
 */
863
void task_get_full_name(int type, int subtype, char *name) {
864
865
866

#ifdef SWIFT_DEBUG_CHECKS
  /* Check input */
867
  if (type >= task_type_count) error("Unknown task type %i", type);
868

869
  if (subtype >= task_subtype_count)
870
871
872
873
874
875
876
877
878
879
    error("Unknown task subtype %i with type %s", subtype, taskID_names[type]);
#endif

  /* Full task name */
  if (subtype == task_subtype_none)
    sprintf(name, "%s", taskID_names[type]);
  else
    sprintf(name, "%s_%s", taskID_names[type], subtaskID_names[subtype]);
}

880
881
882
883
884
885
886
887
888
#ifdef WITH_MPI
/**
 * @brief Create global communicators for each of the subtasks.
 */
void task_create_mpi_comms(void) {
  for (int i = 0; i < task_subtype_count; i++) {
    MPI_Comm_dup(MPI_COMM_WORLD, &subtaskMPI_comms[i]);
  }
}
889
890
891
892
893
894
895
896
/**
 * @brief Create global communicators for each of the subtasks.
 */
void task_free_mpi_comms(void) {
  for (int i = 0; i < task_subtype_count; i++) {
    MPI_Comm_free(&subtaskMPI_comms[i]);
  }
}
897
#endif
898
899

/**
900
901
 * @brief dump all the tasks of all the known engines into a file for
 * postprocessing.
902
903
904
905
906
907
908
909
910
911
912
913
914
915
 *
 * Dumps the information to a file "thread_info-stepn.dat" where n is the
 * given step value, or "thread_info_MPI-stepn.dat", if we are running
 * under MPI. Note if running under MPIU all the ranks are dumped into this
 * one file, which has an additional field to identify the rank.
 *
 * @param e the #engine
 * @param step the current step.
 */
void task_dump_all(struct engine *e, int step) {

#ifdef SWIFT_DEBUG_TASKS

  /* Need this to convert ticks to seconds. */
916
  const unsigned long long cpufreq = clocks_get_cpufreq();
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943

#ifdef WITH_MPI
  /* Make sure output file is empty, only on one rank. */
  char dumpfile[35];
  snprintf(dumpfile, sizeof(dumpfile), "thread_info_MPI-step%d.dat", step);
  FILE *file_thread;
  if (engine_rank == 0) {
    file_thread = fopen(dumpfile, "w");
    fclose(file_thread);
  }
  MPI_Barrier(MPI_COMM_WORLD);

  for (int i = 0; i < e->nr_nodes; i++) {

    /* Rank 0 decides the index of the writing node, this happens
     * one-by-one. */
    int kk = i;
    MPI_Bcast(&kk, 1, MPI_INT, 0, MPI_COMM_WORLD);

    if (i == engine_rank) {

      /* Open file and position at end. */
      file_thread = fopen(dumpfile, "a");

      /* Add some information to help with the plots and conversion of ticks to
       * seconds. */
      fprintf(file_thread, " %03d 0 0 0 0 %lld %lld %lld %lld %lld 0 0 %lld\n",
Josh Borrow's avatar
Josh Borrow committed
944
945
              engine_rank, (long long int)e->tic_step,
              (long long int)e->toc_step, e->updates, e->g_updates,
946
947
948
              e->s_updates, cpufreq);
      int count = 0;
      for (int l = 0; l < e->sched.nr_tasks; l++) {
949
950
        if (!e->sched.tasks[l].implicit &&
            e->sched.tasks[l].tic > e->tic_step) {
951
952
953
954
          fprintf(
              file_thread, " %03i %i %i %i %i %lli %lli %i %i %i %i %lli %i\n",
              engine_rank, e->sched.tasks[l].rid, e->sched.tasks[l].type,
              e->sched.tasks[l].subtype, (e->sched.tasks[l].cj == NULL),
Josh Borrow's avatar
Josh Borrow committed
955
956
              (long long int)e->sched.tasks[l].tic,
              (long long int)e->sched.tasks[l].toc,
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
              (e->sched.tasks[l].ci != NULL) ? e->sched.tasks[l].ci->hydro.count
                                             : 0,
              (e->sched.tasks[l].cj != NULL) ? e->sched.tasks[l].cj->hydro.count
                                             : 0,
              (e->sched.tasks[l].ci != NULL) ? e->sched.tasks[l].ci->grav.count
                                             : 0,
              (e->sched.tasks[l].cj != NULL) ? e->sched.tasks[l].cj->grav.count
                                             : 0,
              e->sched.tasks[l].flags, e->sched.tasks[l].sid);
        }
        count++;
      }
      fclose(file_thread);
    }

    /* And we wait for all to synchronize. */
    MPI_Barrier(MPI_COMM_WORLD);
  }

#else
  /* Non-MPI, so just a single engine's worth of tasks to dump. */
  char dumpfile[32];
  snprintf(dumpfile, sizeof(dumpfile), "thread_info-step%d.dat", step);
  FILE *file_thread;
  file_thread = fopen(dumpfile, "w");

  /* Add some information to help with the plots and conversion of ticks to
   * seconds. */
  fprintf(file_thread, " %d %d %d %d %lld %lld %lld %lld %lld %d %lld\n", -2,
986
987
          -1, -1, 1, (unsigned long long)e->tic_step,
          (unsigned long long)e->toc_step, e->updates, e->g_updates,
988
989
          e->s_updates, 0, cpufreq);
  for (int l = 0; l < e->sched.nr_tasks; l++) {
990
    if (!e->sched.tasks[l].implicit && e->sched.tasks[l].tic > e->tic_step) {
991
992
993
994
      fprintf(
          file_thread, " %i %i %i %i %lli %lli %i %i %i %i %i\n",
          e->sched.tasks[l].rid, e->sched.tasks[l].type,
          e->sched.tasks[l].subtype, (e->sched.tasks[l].cj == NULL),
Josh Borrow's avatar
Josh Borrow committed
995
996
          (unsigned long long)e->sched.tasks[l].tic,
          (unsigned long long)e->sched.tasks[l].toc,
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
          (e->sched.tasks[l].ci == NULL) ? 0
                                         : e->sched.tasks[l].ci->hydro.count,
          (e->sched.tasks[l].cj == NULL) ? 0
                                         : e->sched.tasks[l].cj->hydro.count,
          (e->sched.tasks[l].ci == NULL) ? 0 : e->sched.tasks[l].ci->grav.count,
          (e->sched.tasks[l].cj == NULL) ? 0 : e->sched.tasks[l].cj->grav.count,
          e->sched.tasks[l].sid);
    }
  }
  fclose(file_thread);
#endif  // WITH_MPI
1008
1009
1010
1011
1012
#endif  // SWIFT_DEBUG_TASKS
}

/**
 * @brief Generate simple statistics about the times used by the tasks of
1013
1014
1015
 *        all the engines and write these into two format, a human readable
 *        version for debugging and one intented for inclusion as the fixed
 *        costs for repartitioning.
1016
 *
1017
1018
1019
 * Note that when running under MPI all the tasks can be summed into this single
 * file. In the fuller, human readable file, the statistics included are the
 * number of task of each type/subtype followed by the minimum, maximum, mean
1020
1021
 * and total time taken and the same numbers for the start of the task,
 * in millisec and then the fixed costs value.
1022
 *
1023
1024
1025
 * If header is set, only the fixed costs value is written into the output
 * file in a format that is suitable for inclusion in SWIFT (as
 * partition_fixed_costs.h).
1026
 *
1027
 * @param dumpfile name of the file for the output.
1028
 * @param e the #engine
1029
1030
1031
 * @param header whether to write a header include file.
 * @param allranks do the statistics over all ranks, if not just the current
 *                 one, only used if header is false.
1032
 */
Peter W. Draper's avatar
Peter W. Draper committed
1033
1034
void task_dump_stats(const char *dumpfile, struct engine *e, int header,
                     int allranks) {
1035
1036
1037

  /* Need arrays for sum, min and max across all types and subtypes. */
  double sum[task_type_count][task_subtype_count];
1038
  double tsum[task_type_count][task_subtype_count];
1039
  double min[task_type_count][task_subtype_count];
1040
  double tmin[task_type_count][task_subtype_count];
1041
  double max[task_type_count][task_subtype_count];
1042
  double tmax[task_type_count][task_subtype_count];
1043
1044
1045
1046
1047
  int count[task_type_count][task_subtype_count];

  for (int j = 0; j < task_type_count; j++) {
    for (int k = 0; k < task_subtype_count; k++) {
      sum[j][k] = 0.0;
1048
      tsum[j][k] = 0.0;
1049
1050
      count[j][k] = 0;
      min[j][k] = DBL_MAX;
1051
      tmin[j][k] = DBL_MAX;
1052
      max[j][k] = 0.0;
1053
      tmax[j][k] = 0.0;
1054
1055
1056
1057
1058
1059
1060
    }
  }

  double total[1] = {0.0};
  for (int l = 0; l < e->sched.nr_tasks; l++) {
    int type = e->sched.tasks[l].type;

1061
1062
    /* Skip implicit tasks, tasks that didn't run this step. */
    if (!e->sched.tasks[l].implicit && e->sched.tasks[l].tic > e->tic_step) {
1063
1064
1065
1066
      int subtype = e->sched.tasks[l].subtype;

      double dt = e->sched.tasks[l].toc - e->sched.tasks[l].tic;
      sum[type][subtype] += dt;
1067
1068
1069

      double tic = (double)e->sched.tasks[l].tic;
      tsum[type][subtype] += tic;
1070
1071
1072
1073
      count[type][subtype] += 1;<