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
75
                                             "limiter_in",
                                             "limiter_out",
76
77
78
79
80
81
82
                                             "send",
                                             "recv",
                                             "grav_long_range",
                                             "grav_mm",
                                             "grav_down_in",
                                             "grav_down",
                                             "grav_mesh",
83
                                             "grav_end_force",
84
85
                                             "cooling",
                                             "star_formation",
86
87
                                             "star_formation_in",
                                             "star_formation_out",
88
                                             "logger",
89
90
                                             "stars_in",
                                             "stars_out",
91
92
                                             "stars_ghost_in",
                                             "stars_ghost",
James Willis's avatar
James Willis committed
93
                                             "stars_ghost_out",
94
                                             "stars_sort",
95
                                             "stars_resort",
96
97
                                             "bh_in",
                                             "bh_out",
98
99
100
                                             "bh_density_ghost",
                                             "bh_swallow_ghost1",
                                             "bh_swallow_ghost2",
101
                                             "bh_swallow_ghost3",
102
                                             "fof_self",
James Willis's avatar
James Willis committed
103
                                             "fof_pair"};
104

105
/* Sub-task type names. */
106
const char *subtaskID_names[task_subtype_count] = {
107
108
109
110
111
112
113
114
    "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"};
115

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

121
122
/**
 * @brief Computes the overlap between the parts array of two given cells.
123
 *
Matthieu Schaller's avatar
Matthieu Schaller committed
124
125
126
 * @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.
127
 */
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
#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;                                                               \
  }
145

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

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

  switch (t->type) {

    case task_type_none:
      return task_action_none;
      break;

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

175
176
177
    case task_type_star_formation:
      return task_action_all;

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

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

191
192
193
194
195
196
197
    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:
198
        case task_subtype_gradient:
199
        case task_subtype_force:
200
        case task_subtype_limiter:
201
202
203
          return task_action_part;
          break;

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

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

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

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

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

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

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

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

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

278
279
280
281
282
#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. */
283
  return task_action_none;
284
285
}

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

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

301
302
  /* First check if any of the two tasks are of a type that don't
     use cells. */
303
304
305
306
307
  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);
308
309
  const int ta_spart =
      (ta_act == task_action_spart || ta_act == task_action_all);
310
311
  const int ta_bpart =
      (ta_act == task_action_bpart || ta_act == task_action_all);
312
313
314
  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);
315
316
  const int tb_spart =
      (tb_act == task_action_spart || tb_act == task_action_all);
317
318
  const int tb_bpart =
      (tb_act == task_action_bpart || tb_act == task_action_all);
319
320
321
322
323
324

  /* 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;
325
326
327
328
    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;
329

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

332
333
334
335
336
337
338
339
340
341
342
343
344
345
    /* 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;
346
347
348
349
    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;
350

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

353
354
355
356
357
358
359
360
    /* 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);
  }
361

Loic Hausammann's avatar
Loic Hausammann committed
362
363
364
365
366
  /* 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;
367
368
369
370
    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
371

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

Loic Hausammann's avatar
Loic Hausammann committed
374
375
376
377
378
379
380
381
382
    /* 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);
  }

383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
  /* 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);
  }

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

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

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

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

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

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

439
    case task_type_drift_gpart:
440
    case task_type_grav_mesh:
441
    case task_type_end_grav_force:
442
443
444
      cell_gunlocktree(ci);
      break;

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

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

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

502
    case task_type_grav_down:
503
      cell_gunlocktree(ci);
504
505
506
      cell_munlocktree(ci);
      break;

507
    case task_type_grav_long_range:
508
      cell_munlocktree(ci);
509
      break;
510

511
512
513
514
515
    case task_type_grav_mm:
      cell_munlocktree(ci);
      cell_munlocktree(cj);
      break;

516
517
518
519
    case task_type_star_formation:
      cell_unlocktree(ci);
      cell_sunlocktree(ci);
      cell_gunlocktree(ci);
520
      break;
521

522
523
524
525
    default:
      break;
  }
}
526
527
528
529
530
531

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

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

542
  switch (type) {
543

544
545
546
    /* Communication task? */
    case task_type_recv:
    case task_type_send:
547
#ifdef WITH_MPI
548
549
550
551
552
      /* 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);
553
554
555
556
        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);
557
      }
Peter W. Draper's avatar
Peter W. Draper committed
558
559

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

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

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

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

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

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

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

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

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

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

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

767
768
    default:
      break;
769
770
771
772
773
  }

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

775
776
777
778
779
780
781
782
783
784
785
/**
 * @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);
}
786

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

799
800
  if (type == task_type_grav_long_range || type == task_type_grav_mm ||
      type == task_type_grav_mesh) {
801
802
803
804
805

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

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

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

866
  if (subtype >= task_subtype_count)
867
868
869
870
871
872
873
874
875
876
    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]);
}

877
878
879
880
881
882
883
884
885
#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]);
  }
}
886
887
888
889
890
891
892
893
/**
 * @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]);
  }
}
894
#endif
895
896

/**
897
898
 * @brief dump all the tasks of all the known engines into a file for
 * postprocessing.
899
900
901
902
903
904
905
906
907
908
909
910
911
912
 *
 * 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. */
913
  const unsigned long long cpufreq = clocks_get_cpufreq();
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940

#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
941
942
              engine_rank, (long long int)e->tic_step,
              (long long int)e->toc_step, e->updates, e->g_updates,
943
944
945
              e->s_updates, cpufreq);
      int count = 0;
      for (int l = 0; l < e->sched.nr_tasks; l++) {
946
947
        if (!e->sched.tasks[l].implicit &&
            e->sched.tasks[l].tic > e->tic_step) {
948
949
950
951
          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
952
953
              (long long int)e->sched.tasks[l].tic,
              (long long int)e->sched.tasks[l].toc,
954
955
956
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
              (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,
983
984
          -1, -1, 1, (unsigned long long)e->tic_step,
          (unsigned long long)e->toc_step, e->updates, e->g_updates,
985
986
          e->s_updates, 0, cpufreq);
  for (int l = 0; l < e->sched.nr_tasks; l++) {
987
    if (!e->sched.tasks[l].implicit && e->sched.tasks[l].tic > e->tic_step) {
988
989
990
991
      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
992
993
          (unsigned long long)e->sched.tasks[l].tic,
          (unsigned long long)e->sched.tasks[l].toc,
994
995
996
997
998
999
1000
1001
1002
1003
1004
          (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
1005
1006
1007
1008
1009
#endif  // SWIFT_DEBUG_TASKS
}

/**
 * @brief Generate simple statistics about the times used by the tasks of
1010
1011
1012
 *        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.
1013
 *
1014
1015
1016
 * 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
1017
1018
 * and total time taken and the same numbers for the start of the task,
 * in millisec and then the fixed costs value.
1019
 *
1020
1021
1022
 * 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).
1023
 *
1024
 * @param dumpfile name of the file for the output.
1025
 * @param e the #engine
1026
1027
1028
 * @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.
1029
 */
Peter W. Draper's avatar
Peter W. Draper committed
1030
1031
void task_dump_stats(const char *dumpfile, struct engine *e, int header,
                     int allranks) {
1032
1033
1034

  /* Need arrays for sum, min and max across all types and subtypes. */
  double sum[task_type_count][task_subtype_count];
1035
  double tsum[task_type_count][task_subtype_count];
1036
  double min[task_type_count][task_subtype_count];
1037
  double tmin[task_type_count][task_subtype_count];
1038
  double max[task_type_count][task_subtype_count];
1039
  double tmax[task_type_count][task_subtype_count];
1040
1041
1042
1043
1044
  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;
1045
      tsum[j][k] = 0.0;
1046
1047
      count[j][k] = 0;
      min[j][k] = DBL_MAX;
1048
      tmin[j][k] = DBL_MAX;
1049
      max[j][k] = 0.0;
1050
      tmax[j][k] = 0.0;
1051
1052
1053
1054
1055
1056
1057
    }
  }

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

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

      double dt = e->sched.tasks[l].toc - e->sched.tasks[l].tic;
      sum[type][subtype] += dt;
1064
1065
1066

      double tic = (double)e->sched.tasks[l].tic;
      tsum[type][subtype] += tic;
1067
1068
1069
1070
      count[type][subtype] += 1;
      if (dt < min[type][subtype]) {
        min[type][subtype] = dt;
      }
1071
1072
1073
      if (tic < tmin[type][subtype]) {
        tmin[type][subtype] = tic;
      }
1074
1075
1076