diff --git a/src/cell.c b/src/cell.c
index 232dc14d76689e956bdad2db7d9c48f30165cb7a..c9579be96991b41bfc85eacf0d04176f5994d24b 100644
--- a/src/cell.c
+++ b/src/cell.c
@@ -1059,6 +1059,38 @@ void cell_clear_limiter_flags(struct cell *c, void *data) {
                   cell_flag_do_hydro_limiter | cell_flag_do_hydro_sub_limiter);
 }
 
+void cell_clear_unskip_flags(struct cell *c) {
+
+#ifdef SWIFT_DEBUG_CHECKS
+
+  if (c->split) {
+    for (int k = 0; k < 8; ++k) {
+      if (c->progeny[k] != NULL) cell_clear_unskip_flags(c->progeny[k]);
+    }
+  }
+
+  cell_clear_flag(
+      c, cell_flag_do_stars_resort | cell_flag_do_stars_drift |
+             cell_flag_do_stars_sub_drift | cell_flag_do_hydro_drift |
+             cell_flag_do_hydro_sub_drift | cell_flag_do_hydro_sync |
+             cell_flag_do_hydro_sub_sync | cell_flag_do_grav_drift |
+             cell_flag_do_grav_sub_drift | cell_flag_do_bh_drift |
+             cell_flag_do_bh_sub_drift | cell_flag_do_sink_drift |
+             cell_flag_do_sink_sub_drift | cell_flag_do_hydro_limiter |
+             cell_flag_do_hydro_sub_limiter | cell_flag_do_hydro_sub_sort |
+             cell_flag_do_stars_sub_sort | cell_flag_do_rt_sub_sort |
+             cell_flag_unskip_self_grav_processed |
+             cell_flag_unskip_pair_grav_processed);
+
+  c->hydro.do_sort = 0;
+  c->stars.do_sort = 0;
+  c->hydro.requires_sorts = 0;
+  c->stars.requires_sorts = 0;
+#else
+  error("Calling debugging code without debugging flag activated.");
+#endif
+}
+
 /**
  * @brief Set the super-cell pointers for all cells in a hierarchy.
  *
diff --git a/src/cell.h b/src/cell.h
index 7e9de839e6b87cb8a7f73f07e706b5985cf8b527..eb2c7adeb1f113b8c6f1c248b10cdd42b791998c 100644
--- a/src/cell.h
+++ b/src/cell.h
@@ -637,6 +637,7 @@ void cell_check_spart_pos(const struct cell *c,
 void cell_check_sort_flags(const struct cell *c);
 void cell_clear_stars_sort_flags(struct cell *c, const int unused_flags);
 void cell_clear_hydro_sort_flags(struct cell *c, const int unused_flags);
+void cell_clear_unskip_flags(struct cell *c);
 int cell_has_tasks(struct cell *c);
 void cell_remove_part(const struct engine *e, struct cell *c, struct part *p,
                       struct xpart *xp);
diff --git a/src/cell_black_holes.h b/src/cell_black_holes.h
index fb8a779b4b7c33b8db6d59b44f6480307e2fb36c..01ed209c9d7d988613f13f6efaa80e647ad1668f 100644
--- a/src/cell_black_holes.h
+++ b/src/cell_black_holes.h
@@ -55,9 +55,9 @@ struct cell_black_holes {
     struct task *density_ghost;
 
     /*! The ghost tasks related to BH swallowing */
-    struct task *swallow_ghost_0;
     struct task *swallow_ghost_1;
     struct task *swallow_ghost_2;
+    struct task *swallow_ghost_3;
 
     /*! Linked list of the tasks computing this cell's BH density. */
     struct link *density;
diff --git a/src/cell_unskip.c b/src/cell_unskip.c
index 406505fd8c4ffd66c7f19b0f3331d0f77b29e5e5..fc89d6346ce8238e6d103b999cd0deab319cbb69 100644
--- a/src/cell_unskip.c
+++ b/src/cell_unskip.c
@@ -1087,15 +1087,17 @@ void cell_activate_subcell_black_holes_tasks(struct cell *ci, struct cell *cj,
 
   /* Otherwise, pair interation */
   else {
-    /* Should we even bother? */
-    if (!cell_is_active_black_holes(ci, e) &&
-        !cell_is_active_black_holes(cj, e))
-      return;
 
     /* Get the orientation of the pair. */
     double shift[3];
     const int sid = space_getsid(s->space, &ci, &cj, shift);
 
+    const int ci_active = cell_is_active_black_holes(ci, e);
+    const int cj_active = cell_is_active_black_holes(cj, e);
+
+    /* Should we even bother? */
+    if (!ci_active && !cj_active) return;
+
     /* recurse? */
     if (cell_can_recurse_in_pair_black_holes_task(ci, cj) &&
         cell_can_recurse_in_pair_black_holes_task(cj, ci)) {
@@ -1110,19 +1112,24 @@ void cell_activate_subcell_black_holes_tasks(struct cell *ci, struct cell *cj,
     }
 
     /* Otherwise, activate the drifts. */
-    else if (cell_is_active_black_holes(ci, e) ||
-             cell_is_active_black_holes(cj, e)) {
+    else if (ci_active || cj_active) {
+
+      /* Note we need to drift *both* BH cells to deal with BH<->BH swallows
+       * But we only need to drift the gas cell if the *other* cell has an
+       * active BH */
 
       /* Activate the drifts if the cells are local. */
       if (ci->nodeID == engine_rank) cell_activate_drift_bpart(ci, s);
-      if (cj->nodeID == engine_rank) cell_activate_drift_part(cj, s);
-      if (cj->nodeID == engine_rank && with_timestep_sync)
+      if (cj->nodeID == engine_rank && ci_active)
+        cell_activate_drift_part(cj, s);
+      if (cj->nodeID == engine_rank && ci_active && with_timestep_sync)
         cell_activate_sync_part(cj, s);
 
       /* Activate the drifts if the cells are local. */
-      if (ci->nodeID == engine_rank) cell_activate_drift_part(ci, s);
+      if (ci->nodeID == engine_rank && cj_active)
+        cell_activate_drift_part(ci, s);
       if (cj->nodeID == engine_rank) cell_activate_drift_bpart(cj, s);
-      if (ci->nodeID == engine_rank && with_timestep_sync)
+      if (ci->nodeID == engine_rank && cj_active && with_timestep_sync)
         cell_activate_sync_part(ci, s);
     }
   } /* Otherwise, pair interation */
@@ -2563,7 +2570,8 @@ int cell_unskip_black_holes_tasks(struct cell *c, struct scheduler *s) {
   const int nodeID = e->nodeID;
   int rebuild = 0;
 
-  if (c->black_holes.drift != NULL && cell_is_active_black_holes(c, e)) {
+  if (c->black_holes.drift != NULL && c->black_holes.count > 0 &&
+      cell_is_active_black_holes(c, e)) {
     cell_activate_drift_bpart(c, s);
   }
 
@@ -2572,8 +2580,11 @@ int cell_unskip_black_holes_tasks(struct cell *c, struct scheduler *s) {
     struct task *t = l->t;
     struct cell *ci = t->ci;
     struct cell *cj = t->cj;
-    const int ci_active = cell_is_active_black_holes(ci, e);
-    const int cj_active = (cj != NULL) ? cell_is_active_black_holes(cj, e) : 0;
+    const int ci_active =
+        ci->black_holes.count > 0 && cell_is_active_black_holes(ci, e);
+    const int cj_active = (cj != NULL) ? (cj->black_holes.count > 0 &&
+                                          cell_is_active_black_holes(cj, e))
+                                       : 0;
 #ifdef WITH_MPI
     const int ci_nodeID = ci->nodeID;
     const int cj_nodeID = (cj != NULL) ? cj->nodeID : -1;
@@ -2588,7 +2599,7 @@ int cell_unskip_black_holes_tasks(struct cell *c, struct scheduler *s) {
 
       scheduler_activate(s, t);
 
-      /* Activate the drifts */
+      /* Activate the drifts & sync */
       if (t->type == task_type_self) {
         cell_activate_drift_part(ci, s);
         cell_activate_drift_bpart(ci, s);
@@ -2598,12 +2609,20 @@ int cell_unskip_black_holes_tasks(struct cell *c, struct scheduler *s) {
       /* Activate the drifts */
       else if (t->type == task_type_pair) {
 
-        /* Activate the drift tasks. */
+        /* Activate the drift & sync tasks.
+         * Note we need to drift *both* BH cells to deal with BH<->BH swallows
+         * But we only need to drift the gas cell if the *other* cell has an
+         * active BH */
         if (ci_nodeID == nodeID) cell_activate_drift_bpart(ci, s);
-        if (ci_nodeID == nodeID) cell_activate_drift_part(ci, s);
+        if (ci_nodeID == nodeID && cj_active) cell_activate_drift_part(ci, s);
 
-        if (cj_nodeID == nodeID) cell_activate_drift_part(cj, s);
+        if (cj_nodeID == nodeID && ci_active) cell_activate_drift_part(cj, s);
         if (cj_nodeID == nodeID) cell_activate_drift_bpart(cj, s);
+
+        if (ci_nodeID == nodeID && cj_active && with_timestep_sync)
+          cell_activate_sync_part(ci, s);
+        if (cj_nodeID == nodeID && ci_active && with_timestep_sync)
+          cell_activate_sync_part(cj, s);
       }
 
       /* Store current values of dx_max and h_max. */
@@ -2616,89 +2635,121 @@ int cell_unskip_black_holes_tasks(struct cell *c, struct scheduler *s) {
       else if (t->type == task_type_sub_pair) {
         cell_activate_subcell_black_holes_tasks(ci, cj, s, with_timestep_sync);
       }
+
+      if (t->type == task_type_pair || t->type == task_type_sub_pair) {
+        /* Activate bh_in for each cell that is part of
+         * a pair task as to not miss any dependencies */
+        if (ci_nodeID == nodeID)
+          scheduler_activate(s, ci->hydro.super->black_holes.black_holes_in);
+        if (cj_nodeID == nodeID)
+          scheduler_activate(s, cj->hydro.super->black_holes.black_holes_in);
+      }
     }
 
     /* Only interested in pair interactions as of here. */
     if (t->type == task_type_pair || t->type == task_type_sub_pair) {
 
-      /* Activate bh_in for each cell that is part of
-       * a pair task as to not miss any dependencies */
-      if (ci_nodeID == nodeID)
-        scheduler_activate(s, ci->hydro.super->black_holes.black_holes_in);
-      if (cj_nodeID == nodeID)
-        scheduler_activate(s, cj->hydro.super->black_holes.black_holes_in);
-
       /* Check whether there was too much particle motion, i.e. the
          cell neighbour conditions were violated. */
       if (cell_need_rebuild_for_black_holes_pair(ci, cj)) rebuild = 1;
       if (cell_need_rebuild_for_black_holes_pair(cj, ci)) rebuild = 1;
 
-      scheduler_activate(s, ci->hydro.super->black_holes.swallow_ghost_0);
-      scheduler_activate(s, cj->hydro.super->black_holes.swallow_ghost_0);
+      if (ci->hydro.super->black_holes.count > 0 && ci_active)
+        scheduler_activate(s, ci->hydro.super->black_holes.swallow_ghost_1);
+      if (cj->hydro.super->black_holes.count > 0 && cj_active)
+        scheduler_activate(s, cj->hydro.super->black_holes.swallow_ghost_1);
 
 #ifdef WITH_MPI
       /* Activate the send/recv tasks. */
       if (ci_nodeID != nodeID) {
 
-        /* Receive the foreign parts to compute BH accretion rates and do the
-         * swallowing */
-        scheduler_activate_recv(s, ci->mpi.recv, task_subtype_rho);
-        scheduler_activate_recv(s, ci->mpi.recv, task_subtype_part_swallow);
+        if (ci_active || cj_active) {
+          /* We must exchange the foreign BHs no matter the activity status */
+          scheduler_activate_recv(s, ci->mpi.recv, task_subtype_bpart_rho);
+          scheduler_activate_send(s, cj->mpi.send, task_subtype_bpart_rho,
+                                  ci_nodeID);
 
-        /* Send the local BHs to tag the particles to swallow and to do feedback
-         */
-        scheduler_activate_send(s, cj->mpi.send, task_subtype_bpart_rho,
-                                ci_nodeID);
-        scheduler_activate_send(s, cj->mpi.send, task_subtype_bpart_feedback,
-                                ci_nodeID);
+          /* Drift before you send */
+          if (cj->black_holes.count > 0) cell_activate_drift_bpart(cj, s);
+        }
 
-        /* Drift before you send */
-        cell_activate_drift_bpart(cj, s);
+        if (cj_active) {
 
-        /* Receive the foreign BHs to tag particles to swallow and for feedback
-         */
-        scheduler_activate_recv(s, ci->mpi.recv, task_subtype_bpart_rho);
-        scheduler_activate_recv(s, ci->mpi.recv, task_subtype_bpart_feedback);
+          /* Receive the foreign parts to compute BH accretion rates and do the
+           * swallowing */
+          scheduler_activate_recv(s, ci->mpi.recv, task_subtype_rho);
+          scheduler_activate_recv(s, ci->mpi.recv, task_subtype_part_swallow);
+          scheduler_activate_recv(s, ci->mpi.recv, task_subtype_bpart_merger);
 
-        /* Send the local part information */
-        scheduler_activate_send(s, cj->mpi.send, task_subtype_rho, ci_nodeID);
-        scheduler_activate_send(s, cj->mpi.send, task_subtype_part_swallow,
-                                ci_nodeID);
+          /* Send the local BHs to do feedback */
+          scheduler_activate_send(s, cj->mpi.send, task_subtype_bpart_feedback,
+                                  ci_nodeID);
 
-        /* Drift the cell which will be sent; note that not all sent
-           particles will be drifted, only those that are needed. */
-        cell_activate_drift_part(cj, s);
+          /* Drift before you send */
+          cell_activate_drift_bpart(cj, s);
+        }
+
+        if (ci_active) {
+
+          /* Receive the foreign BHs for feedback */
+          scheduler_activate_recv(s, ci->mpi.recv, task_subtype_bpart_feedback);
+
+          /* Send the local part information */
+          scheduler_activate_send(s, cj->mpi.send, task_subtype_rho, ci_nodeID);
+          scheduler_activate_send(s, cj->mpi.send, task_subtype_part_swallow,
+                                  ci_nodeID);
+          scheduler_activate_send(s, cj->mpi.send, task_subtype_bpart_merger,
+                                  ci_nodeID);
+
+          /* Drift the cell which will be sent; note that not all sent
+             particles will be drifted, only those that are needed. */
+          if (cj->hydro.count > 0) cell_activate_drift_part(cj, s);
+        }
 
       } else if (cj_nodeID != nodeID) {
 
-        /* Receive the foreign parts to compute BH accretion rates and do the
-         * swallowing */
-        scheduler_activate_recv(s, cj->mpi.recv, task_subtype_rho);
-        scheduler_activate_recv(s, cj->mpi.recv, task_subtype_part_swallow);
+        if (ci_active || cj_active) {
+          /* We must exchange the foreign BHs no matter the activity status */
+          scheduler_activate_recv(s, cj->mpi.recv, task_subtype_bpart_rho);
+          scheduler_activate_send(s, ci->mpi.send, task_subtype_bpart_rho,
+                                  cj_nodeID);
 
-        /* Send the local BHs to tag the particles to swallow and to do feedback
-         */
-        scheduler_activate_send(s, ci->mpi.send, task_subtype_bpart_rho,
-                                cj_nodeID);
-        scheduler_activate_send(s, ci->mpi.send, task_subtype_bpart_feedback,
-                                cj_nodeID);
+          /* Drift before you send */
+          if (ci->black_holes.count > 0) cell_activate_drift_bpart(ci, s);
+        }
 
-        /* Drift before you send */
-        cell_activate_drift_bpart(ci, s);
+        if (ci_active) {
 
-        /* Receive the foreign BHs to tag particles to swallow and for feedback
-         */
-        scheduler_activate_recv(s, cj->mpi.recv, task_subtype_bpart_rho);
-        scheduler_activate_recv(s, cj->mpi.recv, task_subtype_bpart_feedback);
+          /* Receive the foreign parts to compute BH accretion rates and do the
+           * swallowing */
+          scheduler_activate_recv(s, cj->mpi.recv, task_subtype_rho);
+          scheduler_activate_recv(s, cj->mpi.recv, task_subtype_part_swallow);
+          scheduler_activate_recv(s, cj->mpi.recv, task_subtype_bpart_merger);
 
-        /* Send the local part information */
-        scheduler_activate_send(s, ci->mpi.send, task_subtype_rho, cj_nodeID);
-        scheduler_activate_send(s, ci->mpi.send, task_subtype_part_swallow,
-                                cj_nodeID);
+          /* Send the local BHs to do feedback */
+          scheduler_activate_send(s, ci->mpi.send, task_subtype_bpart_feedback,
+                                  cj_nodeID);
 
-        /* Drift the cell which will be sent; note that not all sent
-           particles will be drifted, only those that are needed. */
-        cell_activate_drift_part(ci, s);
+          /* Drift before you send */
+          cell_activate_drift_bpart(ci, s);
+        }
+
+        if (cj_active) {
+
+          /* Receive the foreign BHs for feedback */
+          scheduler_activate_recv(s, cj->mpi.recv, task_subtype_bpart_feedback);
+
+          /* Send the local part information */
+          scheduler_activate_send(s, ci->mpi.send, task_subtype_rho, cj_nodeID);
+          scheduler_activate_send(s, ci->mpi.send, task_subtype_part_swallow,
+                                  cj_nodeID);
+          scheduler_activate_send(s, ci->mpi.send, task_subtype_bpart_merger,
+                                  cj_nodeID);
+
+          /* Drift the cell which will be sent; note that not all sent
+             particles will be drifted, only those that are needed. */
+          if (ci->hydro.count > 0) cell_activate_drift_part(ci, s);
+        }
       }
 #endif
     }
@@ -2806,22 +2857,27 @@ int cell_unskip_black_holes_tasks(struct cell *c, struct scheduler *s) {
   }
 
   /* Unskip all the other task types. */
-  if (c->nodeID == nodeID && cell_is_active_black_holes(c, e)) {
-
-    /* If the cell doesn't have any pair/sub_pair type tasks,
-     * then we haven't unskipped all the implicit tasks yet. */
+  if (cell_is_active_black_holes(c, e)) {
     if (c->black_holes.density_ghost != NULL)
       scheduler_activate(s, c->black_holes.density_ghost);
-    if (c->black_holes.swallow_ghost_0 != NULL)
-      scheduler_activate(s, c->black_holes.swallow_ghost_0);
     if (c->black_holes.swallow_ghost_1 != NULL)
       scheduler_activate(s, c->black_holes.swallow_ghost_1);
     if (c->black_holes.swallow_ghost_2 != NULL)
       scheduler_activate(s, c->black_holes.swallow_ghost_2);
+    if (c->black_holes.swallow_ghost_3 != NULL)
+      scheduler_activate(s, c->black_holes.swallow_ghost_3);
+  }
+  if (c->nodeID == nodeID && cell_is_active_black_holes(c, e)) {
     if (c->black_holes.black_holes_in != NULL)
       scheduler_activate(s, c->black_holes.black_holes_in);
     if (c->black_holes.black_holes_out != NULL)
       scheduler_activate(s, c->black_holes.black_holes_out);
+  }
+  if (c->nodeID == nodeID && c->black_holes.count > 0 &&
+      cell_is_active_black_holes(c, e)) {
+
+    /* If the cell doesn't have any pair/sub_pair type tasks,
+     * then we haven't unskipped all the implicit tasks yet. */
     if (c->kick1 != NULL) scheduler_activate(s, c->kick1);
     if (c->kick2 != NULL) scheduler_activate(s, c->kick2);
     if (c->timestep != NULL) scheduler_activate(s, c->timestep);
diff --git a/src/engine.c b/src/engine.c
index 83d2c3cbb0df4fed49621e38cc9dda59f675e44b..35d60b0281a8c5b2ab5b26f670bce614a0f5bf81 100644
--- a/src/engine.c
+++ b/src/engine.c
@@ -144,6 +144,10 @@ int engine_rank;
 /** The current step of the engine as a global variable (for messages). */
 int engine_current_step;
 
+#ifdef SWIFT_DEBUG_CHECKS
+extern int activate_by_unskip;
+#endif
+
 /**
  * @brief Link a density/force task to a cell.
  *
@@ -1363,8 +1367,48 @@ void engine_rebuild(struct engine *e, const int repartitioned,
 #endif
 
   /* Run through the tasks and mark as skip or not. */
-  if (engine_marktasks(e))
-    error("engine_marktasks failed after space_rebuild.");
+#ifdef SWIFT_DEBUG_CHECKS
+  activate_by_unskip = 1;
+#endif
+  engine_unskip(e);
+  if (e->forcerebuild) error("engine_unskip faled after a rebuild!");
+
+#ifdef SWIFT_DEBUG_CHECKS
+
+  /* Reset all the tasks */
+  for (int i = 0; i < e->sched.nr_tasks; ++i) {
+    e->sched.tasks[i].skip = 1;
+  }
+  for (int i = 0; i < e->sched.active_count; ++i) {
+    e->sched.tid_active[i] = -1;
+  }
+  e->sched.active_count = 0;
+  for (int i = 0; i < e->s->nr_cells; ++i) {
+    cell_clear_unskip_flags(&e->s->cells_top[i]);
+  }
+
+  /* Now run the (legacy) marktasks */
+  activate_by_unskip = 0;
+  engine_marktasks(e);
+
+  /* Verify that the two task activation procedures match */
+  for (int i = 0; i < e->sched.nr_tasks; ++i) {
+    struct task *t = &e->sched.tasks[i];
+
+    if (t->activated_by_unskip && !t->activated_by_marktask) {
+      error("Task %s/%s activated by unskip and not by marktask!",
+            taskID_names[t->type], subtaskID_names[t->subtype]);
+    }
+
+    if (!t->activated_by_unskip && t->activated_by_marktask) {
+      error("Task %s/%s activated by marktask and not by unskip!",
+            taskID_names[t->type], subtaskID_names[t->subtype]);
+    }
+
+    t->activated_by_marktask = 0;
+    t->activated_by_unskip = 0;
+  }
+#endif
 
   /* Print the status of the system */
   if (e->verbose) engine_print_task_counts(e);
diff --git a/src/engine_maketasks.c b/src/engine_maketasks.c
index 5b5873ace0bdd983e13cfb731415cd82adb1132f..53468bbf5efe40bc617f250f3eb810a483264008 100644
--- a/src/engine_maketasks.c
+++ b/src/engine_maketasks.c
@@ -505,24 +505,24 @@ void engine_addtasks_send_black_holes(struct engine *e, struct cell *ci,
       scheduler_addunlock(s, t_feedback,
                           ci->hydro.super->black_holes.black_holes_out);
 
-      scheduler_addunlock(s, ci->hydro.super->black_holes.swallow_ghost_2,
+      scheduler_addunlock(s, ci->hydro.super->black_holes.swallow_ghost_3,
                           t_feedback);
 
       /* Ghost before you send */
       scheduler_addunlock(s, ci->hydro.super->black_holes.drift, t_rho);
       scheduler_addunlock(s, ci->hydro.super->black_holes.density_ghost, t_rho);
       scheduler_addunlock(s, t_rho,
-                          ci->hydro.super->black_holes.swallow_ghost_0);
+                          ci->hydro.super->black_holes.swallow_ghost_1);
 
-      scheduler_addunlock(s, ci->hydro.super->black_holes.swallow_ghost_0,
+      scheduler_addunlock(s, ci->hydro.super->black_holes.swallow_ghost_1,
                           t_bh_merger);
       scheduler_addunlock(s, t_bh_merger,
-                          ci->hydro.super->black_holes.swallow_ghost_2);
+                          ci->hydro.super->black_holes.swallow_ghost_3);
 
-      scheduler_addunlock(s, ci->hydro.super->black_holes.swallow_ghost_0,
+      scheduler_addunlock(s, ci->hydro.super->black_holes.swallow_ghost_1,
                           t_gas_swallow);
       scheduler_addunlock(s, t_gas_swallow,
-                          ci->hydro.super->black_holes.swallow_ghost_1);
+                          ci->hydro.super->black_holes.swallow_ghost_2);
     }
 
     engine_addlink(e, &ci->mpi.send, t_rho);
@@ -1559,7 +1559,7 @@ void engine_make_hierarchical_tasks_hydro(struct engine *e, struct cell *c,
     }
 
     if (with_black_holes) {
-      c->black_holes.swallow_ghost_0 =
+      c->black_holes.swallow_ghost_1 =
           scheduler_addtask(s, task_type_bh_swallow_ghost1, task_subtype_none,
                             0, /* implicit =*/1, c, NULL);
     }
@@ -1794,11 +1794,11 @@ void engine_make_hierarchical_tasks_hydro(struct engine *e, struct cell *c,
         c->black_holes.density_ghost = scheduler_addtask(
             s, task_type_bh_density_ghost, task_subtype_none, 0, 0, c, NULL);
 
-        c->black_holes.swallow_ghost_1 =
+        c->black_holes.swallow_ghost_2 =
             scheduler_addtask(s, task_type_bh_swallow_ghost2, task_subtype_none,
                               0, /* implicit =*/1, c, NULL);
 
-        c->black_holes.swallow_ghost_2 = scheduler_addtask(
+        c->black_holes.swallow_ghost_3 = scheduler_addtask(
             s, task_type_bh_swallow_ghost3, task_subtype_none, 0, 0, c, NULL);
 
 #ifdef WITH_CSDS
@@ -1820,7 +1820,7 @@ void engine_make_hierarchical_tasks_hydro(struct engine *e, struct cell *c,
         /* Make sure we don't start swallowing gas particles before the stars
            have converged on their smoothing lengths. */
         scheduler_addunlock(s, c->stars.density_ghost,
-                            c->black_holes.swallow_ghost_0);
+                            c->black_holes.swallow_ghost_1);
       }
     }
   } else { /* We are above the super-cell so need to go deeper */
@@ -2679,19 +2679,19 @@ void engine_make_extra_hydroloop_tasks_mapper(void *map_data, int num_elements,
         scheduler_addunlock(sched, ci->hydro.super->black_holes.density_ghost,
                             t_bh_swallow);
         scheduler_addunlock(sched, t_bh_swallow,
-                            ci->hydro.super->black_holes.swallow_ghost_0);
+                            ci->hydro.super->black_holes.swallow_ghost_1);
 
-        scheduler_addunlock(sched, ci->hydro.super->black_holes.swallow_ghost_0,
+        scheduler_addunlock(sched, ci->hydro.super->black_holes.swallow_ghost_1,
                             t_do_gas_swallow);
         scheduler_addunlock(sched, t_do_gas_swallow,
-                            ci->hydro.super->black_holes.swallow_ghost_1);
+                            ci->hydro.super->black_holes.swallow_ghost_2);
 
-        scheduler_addunlock(sched, ci->hydro.super->black_holes.swallow_ghost_1,
+        scheduler_addunlock(sched, ci->hydro.super->black_holes.swallow_ghost_2,
                             t_do_bh_swallow);
         scheduler_addunlock(sched, t_do_bh_swallow,
-                            ci->hydro.super->black_holes.swallow_ghost_2);
+                            ci->hydro.super->black_holes.swallow_ghost_3);
 
-        scheduler_addunlock(sched, ci->hydro.super->black_holes.swallow_ghost_2,
+        scheduler_addunlock(sched, ci->hydro.super->black_holes.swallow_ghost_3,
                             t_bh_feedback);
         scheduler_addunlock(sched, t_bh_feedback,
                             ci->hydro.super->black_holes.black_holes_out);
@@ -3037,22 +3037,22 @@ void engine_make_extra_hydroloop_tasks_mapper(void *map_data, int num_elements,
           scheduler_addunlock(sched, ci->hydro.super->black_holes.density_ghost,
                               t_bh_swallow);
           scheduler_addunlock(sched, t_bh_swallow,
-                              ci->hydro.super->black_holes.swallow_ghost_0);
+                              ci->hydro.super->black_holes.swallow_ghost_1);
 
           scheduler_addunlock(sched,
-                              ci->hydro.super->black_holes.swallow_ghost_0,
+                              ci->hydro.super->black_holes.swallow_ghost_1,
                               t_do_gas_swallow);
           scheduler_addunlock(sched, t_do_gas_swallow,
-                              ci->hydro.super->black_holes.swallow_ghost_1);
+                              ci->hydro.super->black_holes.swallow_ghost_2);
 
           scheduler_addunlock(sched,
-                              ci->hydro.super->black_holes.swallow_ghost_1,
+                              ci->hydro.super->black_holes.swallow_ghost_2,
                               t_do_bh_swallow);
           scheduler_addunlock(sched, t_do_bh_swallow,
-                              ci->hydro.super->black_holes.swallow_ghost_2);
+                              ci->hydro.super->black_holes.swallow_ghost_3);
 
           scheduler_addunlock(sched,
-                              ci->hydro.super->black_holes.swallow_ghost_2,
+                              ci->hydro.super->black_holes.swallow_ghost_3,
                               t_bh_feedback);
           scheduler_addunlock(sched, t_bh_feedback,
                               ci->hydro.super->black_holes.black_holes_out);
@@ -3098,7 +3098,7 @@ void engine_make_extra_hydroloop_tasks_mapper(void *map_data, int num_elements,
 
         if (with_black_holes && (bcount_i > 0 || bcount_j > 0)) {
           scheduler_addunlock(sched, t_bh_swallow,
-                              ci->hydro.super->black_holes.swallow_ghost_0);
+                              ci->hydro.super->black_holes.swallow_ghost_1);
         }
       }
 
@@ -3196,22 +3196,22 @@ void engine_make_extra_hydroloop_tasks_mapper(void *map_data, int num_elements,
                                 cj->hydro.super->black_holes.density_ghost,
                                 t_bh_swallow);
             scheduler_addunlock(sched, t_bh_swallow,
-                                cj->hydro.super->black_holes.swallow_ghost_0);
+                                cj->hydro.super->black_holes.swallow_ghost_1);
 
             scheduler_addunlock(sched,
-                                cj->hydro.super->black_holes.swallow_ghost_0,
+                                cj->hydro.super->black_holes.swallow_ghost_1,
                                 t_do_gas_swallow);
             scheduler_addunlock(sched, t_do_gas_swallow,
-                                cj->hydro.super->black_holes.swallow_ghost_1);
+                                cj->hydro.super->black_holes.swallow_ghost_2);
 
             scheduler_addunlock(sched,
-                                cj->hydro.super->black_holes.swallow_ghost_1,
+                                cj->hydro.super->black_holes.swallow_ghost_2,
                                 t_do_bh_swallow);
             scheduler_addunlock(sched, t_do_bh_swallow,
-                                cj->hydro.super->black_holes.swallow_ghost_2);
+                                cj->hydro.super->black_holes.swallow_ghost_3);
 
             scheduler_addunlock(sched,
-                                cj->hydro.super->black_holes.swallow_ghost_2,
+                                cj->hydro.super->black_holes.swallow_ghost_3,
                                 t_bh_feedback);
             scheduler_addunlock(sched, t_bh_feedback,
                                 cj->hydro.super->black_holes.black_holes_out);
@@ -3266,7 +3266,7 @@ void engine_make_extra_hydroloop_tasks_mapper(void *map_data, int num_elements,
         if (with_black_holes && (bcount_i > 0 || bcount_j > 0)) {
 
           scheduler_addunlock(sched, t_bh_swallow,
-                              cj->hydro.super->black_holes.swallow_ghost_0);
+                              cj->hydro.super->black_holes.swallow_ghost_1);
         }
       }
     }
@@ -3496,19 +3496,19 @@ void engine_make_extra_hydroloop_tasks_mapper(void *map_data, int num_elements,
         scheduler_addunlock(sched, ci->hydro.super->black_holes.density_ghost,
                             t_bh_swallow);
         scheduler_addunlock(sched, t_bh_swallow,
-                            ci->hydro.super->black_holes.swallow_ghost_0);
+                            ci->hydro.super->black_holes.swallow_ghost_1);
 
-        scheduler_addunlock(sched, ci->hydro.super->black_holes.swallow_ghost_0,
+        scheduler_addunlock(sched, ci->hydro.super->black_holes.swallow_ghost_1,
                             t_do_gas_swallow);
         scheduler_addunlock(sched, t_do_gas_swallow,
-                            ci->hydro.super->black_holes.swallow_ghost_1);
+                            ci->hydro.super->black_holes.swallow_ghost_2);
 
-        scheduler_addunlock(sched, ci->hydro.super->black_holes.swallow_ghost_1,
+        scheduler_addunlock(sched, ci->hydro.super->black_holes.swallow_ghost_2,
                             t_do_bh_swallow);
         scheduler_addunlock(sched, t_do_bh_swallow,
-                            ci->hydro.super->black_holes.swallow_ghost_2);
+                            ci->hydro.super->black_holes.swallow_ghost_3);
 
-        scheduler_addunlock(sched, ci->hydro.super->black_holes.swallow_ghost_2,
+        scheduler_addunlock(sched, ci->hydro.super->black_holes.swallow_ghost_3,
                             t_bh_feedback);
         scheduler_addunlock(sched, t_bh_feedback,
                             ci->hydro.super->black_holes.black_holes_out);
@@ -3866,22 +3866,22 @@ void engine_make_extra_hydroloop_tasks_mapper(void *map_data, int num_elements,
           scheduler_addunlock(sched, ci->hydro.super->black_holes.density_ghost,
                               t_bh_swallow);
           scheduler_addunlock(sched, t_bh_swallow,
-                              ci->hydro.super->black_holes.swallow_ghost_0);
+                              ci->hydro.super->black_holes.swallow_ghost_1);
 
           scheduler_addunlock(sched,
-                              ci->hydro.super->black_holes.swallow_ghost_0,
+                              ci->hydro.super->black_holes.swallow_ghost_1,
                               t_do_gas_swallow);
           scheduler_addunlock(sched, t_do_gas_swallow,
-                              ci->hydro.super->black_holes.swallow_ghost_1);
+                              ci->hydro.super->black_holes.swallow_ghost_2);
 
           scheduler_addunlock(sched,
-                              ci->hydro.super->black_holes.swallow_ghost_1,
+                              ci->hydro.super->black_holes.swallow_ghost_2,
                               t_do_bh_swallow);
           scheduler_addunlock(sched, t_do_bh_swallow,
-                              ci->hydro.super->black_holes.swallow_ghost_2);
+                              ci->hydro.super->black_holes.swallow_ghost_3);
 
           scheduler_addunlock(sched,
-                              ci->hydro.super->black_holes.swallow_ghost_2,
+                              ci->hydro.super->black_holes.swallow_ghost_3,
                               t_bh_feedback);
           scheduler_addunlock(sched, t_bh_feedback,
                               ci->hydro.super->black_holes.black_holes_out);
@@ -3927,7 +3927,7 @@ void engine_make_extra_hydroloop_tasks_mapper(void *map_data, int num_elements,
         if (with_black_holes && (bcount_i > 0 || bcount_j > 0)) {
 
           scheduler_addunlock(sched, t_bh_swallow,
-                              ci->hydro.super->black_holes.swallow_ghost_0);
+                              ci->hydro.super->black_holes.swallow_ghost_1);
         }
       }
 
@@ -4024,22 +4024,22 @@ void engine_make_extra_hydroloop_tasks_mapper(void *map_data, int num_elements,
                                 cj->hydro.super->black_holes.density_ghost,
                                 t_bh_swallow);
             scheduler_addunlock(sched, t_bh_swallow,
-                                cj->hydro.super->black_holes.swallow_ghost_0);
+                                cj->hydro.super->black_holes.swallow_ghost_1);
 
             scheduler_addunlock(sched,
-                                cj->hydro.super->black_holes.swallow_ghost_0,
+                                cj->hydro.super->black_holes.swallow_ghost_1,
                                 t_do_gas_swallow);
             scheduler_addunlock(sched, t_do_gas_swallow,
-                                cj->hydro.super->black_holes.swallow_ghost_1);
+                                cj->hydro.super->black_holes.swallow_ghost_2);
 
             scheduler_addunlock(sched,
-                                cj->hydro.super->black_holes.swallow_ghost_1,
+                                cj->hydro.super->black_holes.swallow_ghost_2,
                                 t_do_bh_swallow);
             scheduler_addunlock(sched, t_do_bh_swallow,
-                                cj->hydro.super->black_holes.swallow_ghost_2);
+                                cj->hydro.super->black_holes.swallow_ghost_3);
 
             scheduler_addunlock(sched,
-                                cj->hydro.super->black_holes.swallow_ghost_2,
+                                cj->hydro.super->black_holes.swallow_ghost_3,
                                 t_bh_feedback);
             scheduler_addunlock(sched, t_bh_feedback,
                                 cj->hydro.super->black_holes.black_holes_out);
@@ -4091,7 +4091,7 @@ void engine_make_extra_hydroloop_tasks_mapper(void *map_data, int num_elements,
 
         if (with_black_holes && (bcount_i > 0 || bcount_j > 0)) {
           scheduler_addunlock(sched, t_bh_swallow,
-                              cj->hydro.super->black_holes.swallow_ghost_0);
+                              cj->hydro.super->black_holes.swallow_ghost_1);
         }
       }
     }
diff --git a/src/engine_marktasks.c b/src/engine_marktasks.c
index a5e6ad5b2ef6697938340b5b4c3c833d28876694..576f11e2d8fa3ea6b91eda6a93dc7bfe0774de85 100644
--- a/src/engine_marktasks.c
+++ b/src/engine_marktasks.c
@@ -97,7 +97,8 @@ void engine_marktasks_mapper(void *map_data, int num_elements,
 
       const int ci_active_hydro = cell_is_active_hydro(ci, e);
       const int ci_active_gravity = cell_is_active_gravity(ci, e);
-      const int ci_active_black_holes = cell_is_active_black_holes(ci, e);
+      const int ci_active_black_holes =
+          ci->black_holes.count > 0 && cell_is_active_black_holes(ci, e);
       const int ci_active_sinks =
           cell_is_active_sinks(ci, e) || ci_active_hydro;
       const int ci_active_stars = cell_need_activating_stars(
@@ -342,7 +343,7 @@ void engine_marktasks_mapper(void *map_data, int num_elements,
                t_subtype == task_subtype_external_grav) {
         if (ci_active_gravity) {
           scheduler_activate(s, t);
-          cell_activate_drift_gpart(t->ci, s);
+          cell_activate_subcell_external_grav_tasks(t->ci, s);
         }
       }
 
@@ -390,8 +391,10 @@ void engine_marktasks_mapper(void *map_data, int num_elements,
       const int ci_active_gravity = cell_is_active_gravity(ci, e);
       const int cj_active_gravity = cell_is_active_gravity(cj, e);
 
-      const int ci_active_black_holes = cell_is_active_black_holes(ci, e);
-      const int cj_active_black_holes = cell_is_active_black_holes(cj, e);
+      const int ci_active_black_holes =
+          ci->black_holes.count > 0 && cell_is_active_black_holes(ci, e);
+      const int cj_active_black_holes =
+          cj->black_holes.count > 0 && cell_is_active_black_holes(cj, e);
 
       const int ci_active_sinks =
           cell_is_active_sinks(ci, e) || ci_active_hydro;
@@ -605,12 +608,15 @@ void engine_marktasks_mapper(void *map_data, int num_elements,
         }
 
         if (t->type == task_type_pair || t->type == task_type_sub_pair) {
-          /* Add stars_out dependencies for each cell that is part of
-           * a pair/sub_pair task as to not miss any dependencies */
-          if (ci_nodeID == nodeID)
-            scheduler_activate(s, ci->hydro.super->stars.stars_out);
-          if (cj_nodeID == nodeID)
-            scheduler_activate(s, cj->hydro.super->stars.stars_out);
+
+          if (ci_active_stars || cj_active_stars) {
+            /* Add stars_out dependencies for each cell that is part of
+             * a pair/sub_pair task as to not miss any dependencies */
+            if (ci_nodeID == nodeID)
+              scheduler_activate(s, ci->hydro.super->stars.stars_out);
+            if (cj_nodeID == nodeID)
+              scheduler_activate(s, cj->hydro.super->stars.stars_out);
+          }
         }
       }
 
@@ -627,12 +633,36 @@ void engine_marktasks_mapper(void *map_data, int num_elements,
 
         /* Set the correct drifting flags */
         if (t_type == task_type_pair && t_subtype == task_subtype_bh_density) {
+
+          /* Note we need to drift *both* BH cells to deal with BH<->BH swallows
+           * But we only need to drift the gas cell if the *other* cell has an
+           * active BH */
           if (ci_nodeID == nodeID) cell_activate_drift_bpart(ci, s);
-          if (ci_nodeID == nodeID) cell_activate_drift_part(ci, s);
+          if (ci_nodeID == nodeID && cj_active_black_holes)
+            cell_activate_drift_part(ci, s);
 
-          if (cj_nodeID == nodeID) cell_activate_drift_part(cj, s);
+          if (cj_nodeID == nodeID && ci_active_black_holes)
+            cell_activate_drift_part(cj, s);
           if (cj_nodeID == nodeID) cell_activate_drift_bpart(cj, s);
 
+          if (ci_nodeID == nodeID && cj_active_black_holes &&
+              with_timestep_sync)
+            cell_activate_sync_part(ci, s);
+          if (cj_nodeID == nodeID && ci_active_black_holes &&
+              with_timestep_sync)
+            cell_activate_sync_part(cj, s);
+        }
+
+        /* Store current values of dx_max and h_max. */
+        else if (t_type == task_type_sub_pair &&
+                 t_subtype == task_subtype_bh_density) {
+          cell_activate_subcell_black_holes_tasks(ci, cj, s,
+                                                  with_timestep_sync);
+        }
+
+        if ((t_type == task_type_pair || t_type == task_type_sub_pair) &&
+            t_subtype == task_subtype_bh_density) {
+
           /* Activate bh_in for each cell that is part of
            * a pair task as to not miss any dependencies */
           if (ci_nodeID == nodeID)
@@ -643,6 +673,7 @@ void engine_marktasks_mapper(void *map_data, int num_elements,
 
         if ((t_type == task_type_pair || t_type == task_type_sub_pair) &&
             t_subtype == task_subtype_bh_feedback) {
+
           /* Add bh_out dependencies for each cell that is part of
            * a pair/sub_pair task as to not miss any dependencies */
           if (ci_nodeID == nodeID)
@@ -650,19 +681,6 @@ void engine_marktasks_mapper(void *map_data, int num_elements,
           if (cj_nodeID == nodeID)
             scheduler_activate(s, cj->hydro.super->black_holes.black_holes_out);
         }
-
-        /* Store current values of dx_max and h_max. */
-        else if (t_type == task_type_sub_pair &&
-                 t_subtype == task_subtype_bh_density) {
-          cell_activate_subcell_black_holes_tasks(ci, cj, s,
-                                                  with_timestep_sync);
-          /* Activate bh_in for each cell that is part of
-           * a sub_pair task as to not miss any dependencies */
-          if (ci_nodeID == nodeID)
-            scheduler_activate(s, ci->hydro.super->black_holes.black_holes_in);
-          if (cj_nodeID == nodeID)
-            scheduler_activate(s, cj->hydro.super->black_holes.black_holes_in);
-        }
       }
 
       /* Gravity */
@@ -871,6 +889,7 @@ void engine_marktasks_mapper(void *map_data, int num_elements,
                                               with_timestep_limiter);
           }
         }
+#endif
       }
 
       /* Pair tasks between inactive local cells and active remote cells. */
@@ -879,6 +898,7 @@ void engine_marktasks_mapper(void *map_data, int num_elements,
           (ci_nodeID == nodeID && cj_nodeID != nodeID && !ci_active_rt &&
            cj_active_rt)) {
 
+#if defined(WITH_MPI) && defined(MPI_SYMMETRIC_FORCE_INTERACTION)
         if (t_subtype == task_subtype_rt_transport) {
 
           scheduler_activate(s, t);
@@ -955,13 +975,14 @@ void engine_marktasks_mapper(void *map_data, int num_elements,
 
           /* Is the foreign cell active and will need stuff from us? */
           if (ci_active_hydro) {
-            struct link *l = scheduler_activate_send(
-                s, cj->mpi.send, task_subtype_xv, ci_nodeID);
+            scheduler_activate_send(s, cj->mpi.send, task_subtype_xv,
+                                    ci_nodeID);
 
             /* Drift the cell which will be sent at the level at which it is
                sent, i.e. drift the cell specified in the send task (l->t)
                itself. */
-            cell_activate_drift_part(l->t->ci, s);
+            cell_activate_drift_part(cj, s);
+            if (with_timestep_limiter) cell_activate_limiter(cj, s);
 
             /* If the local cell is also active, more stuff will be needed. */
             if (cj_active_hydro) {
@@ -1055,13 +1076,14 @@ void engine_marktasks_mapper(void *map_data, int num_elements,
           /* Is the foreign cell active and will need stuff from us? */
           if (cj_active_hydro) {
 
-            struct link *l = scheduler_activate_send(
-                s, ci->mpi.send, task_subtype_xv, cj_nodeID);
+            scheduler_activate_send(s, ci->mpi.send, task_subtype_xv,
+                                    cj_nodeID);
 
             /* Drift the cell which will be sent at the level at which it is
                sent, i.e. drift the cell specified in the send task (l->t)
                itself. */
-            cell_activate_drift_part(l->t->ci, s);
+            cell_activate_drift_part(ci, s);
+            if (with_timestep_limiter) cell_activate_limiter(ci, s);
 
             /* If the local cell is also active, more stuff will be needed. */
             if (ci_active_hydro) {
@@ -1235,78 +1257,106 @@ void engine_marktasks_mapper(void *map_data, int num_elements,
         if (cell_need_rebuild_for_black_holes_pair(ci, cj)) *rebuild_space = 1;
         if (cell_need_rebuild_for_black_holes_pair(cj, ci)) *rebuild_space = 1;
 
-        scheduler_activate(s, ci->hydro.super->black_holes.swallow_ghost_0);
-        scheduler_activate(s, cj->hydro.super->black_holes.swallow_ghost_0);
+        if (ci->hydro.super->black_holes.count > 0 && ci_active_black_holes)
+          scheduler_activate(s, ci->hydro.super->black_holes.swallow_ghost_1);
+        if (cj->hydro.super->black_holes.count > 0 && cj_active_black_holes)
+          scheduler_activate(s, cj->hydro.super->black_holes.swallow_ghost_1);
 
 #ifdef WITH_MPI
         /* Activate the send/recv tasks. */
         if (ci_nodeID != nodeID) {
 
-          /* Receive the foreign parts to compute BH accretion rates and do the
-           * swallowing */
-          scheduler_activate_recv(s, ci->mpi.recv, task_subtype_rho);
-          scheduler_activate_recv(s, ci->mpi.recv, task_subtype_part_swallow);
-          scheduler_activate_recv(s, ci->mpi.recv, task_subtype_bpart_merger);
-
-          /* Send the local BHs to tag the particles to swallow and to do
-           * feedback */
-          scheduler_activate_send(s, cj->mpi.send, task_subtype_bpart_rho,
-                                  ci_nodeID);
-          scheduler_activate_send(s, cj->mpi.send, task_subtype_bpart_feedback,
-                                  ci_nodeID);
-
-          /* Drift before you send */
-          cell_activate_drift_bpart(cj, s);
-
-          /* Receive the foreign BHs to tag particles to swallow and for
-           * feedback */
-          scheduler_activate_recv(s, ci->mpi.recv, task_subtype_bpart_rho);
-          scheduler_activate_recv(s, ci->mpi.recv, task_subtype_bpart_feedback);
-
-          /* Send the local part information */
-          scheduler_activate_send(s, cj->mpi.send, task_subtype_rho, ci_nodeID);
-          scheduler_activate_send(s, cj->mpi.send, task_subtype_part_swallow,
-                                  ci_nodeID);
-          scheduler_activate_send(s, cj->mpi.send, task_subtype_bpart_merger,
-                                  ci_nodeID);
-
-          /* Drift the cell which will be sent; note that not all sent
-             particles will be drifted, only those that are needed. */
-          cell_activate_drift_part(cj, s);
+          if (ci_active_black_holes || cj_active_black_holes) {
+            /* We must exchange the foreign BHs no matter the activity status */
+            scheduler_activate_recv(s, ci->mpi.recv, task_subtype_bpart_rho);
+            scheduler_activate_send(s, cj->mpi.send, task_subtype_bpart_rho,
+                                    ci_nodeID);
+
+            /* Drift before you send */
+            if (cj->black_holes.count > 0) cell_activate_drift_bpart(cj, s);
+          }
+
+          if (cj_active_black_holes) {
+
+            /* Receive the foreign parts to compute BH accretion rates and do
+             * the swallowing */
+            scheduler_activate_recv(s, ci->mpi.recv, task_subtype_rho);
+            scheduler_activate_recv(s, ci->mpi.recv, task_subtype_part_swallow);
+            scheduler_activate_recv(s, ci->mpi.recv, task_subtype_bpart_merger);
+
+            /* Send the local BHs to tag the particles to do feedback */
+            scheduler_activate_send(s, cj->mpi.send,
+                                    task_subtype_bpart_feedback, ci_nodeID);
+
+            /* Drift before you send */
+            cell_activate_drift_bpart(cj, s);
+          }
+
+          if (ci_active_black_holes) {
+
+            /* Receive the foreign BHs for feedback */
+            scheduler_activate_recv(s, ci->mpi.recv,
+                                    task_subtype_bpart_feedback);
+
+            /* Send the local part information */
+            scheduler_activate_send(s, cj->mpi.send, task_subtype_rho,
+                                    ci_nodeID);
+            scheduler_activate_send(s, cj->mpi.send, task_subtype_part_swallow,
+                                    ci_nodeID);
+            scheduler_activate_send(s, cj->mpi.send, task_subtype_bpart_merger,
+                                    ci_nodeID);
+
+            /* Drift the cell which will be sent; note that not all sent
+               particles will be drifted, only those that are needed. */
+            if (cj->hydro.count > 0) cell_activate_drift_part(cj, s);
+          }
 
         } else if (cj_nodeID != nodeID) {
 
-          /* Receive the foreign parts to compute BH accretion rates and do the
-           * swallowing */
-          scheduler_activate_recv(s, cj->mpi.recv, task_subtype_rho);
-          scheduler_activate_recv(s, cj->mpi.recv, task_subtype_part_swallow);
-          scheduler_activate_recv(s, cj->mpi.recv, task_subtype_bpart_merger);
+          if (ci_active_black_holes || cj_active_black_holes) {
+            /* We must exchange the foreign BHs no matter the activity status */
+            scheduler_activate_recv(s, cj->mpi.recv, task_subtype_bpart_rho);
+            scheduler_activate_send(s, ci->mpi.send, task_subtype_bpart_rho,
+                                    cj_nodeID);
 
-          /* Send the local BHs to tag the particles to swallow and to do
-           * feedback */
-          scheduler_activate_send(s, ci->mpi.send, task_subtype_bpart_rho,
-                                  cj_nodeID);
-          scheduler_activate_send(s, ci->mpi.send, task_subtype_bpart_feedback,
-                                  cj_nodeID);
+            /* Drift before you send */
+            if (ci->black_holes.count > 0) cell_activate_drift_bpart(ci, s);
+          }
 
-          /* Drift before you send */
-          cell_activate_drift_bpart(ci, s);
+          if (ci_active_black_holes) {
 
-          /* Receive the foreign BHs to tag particles to swallow and for
-           * feedback */
-          scheduler_activate_recv(s, cj->mpi.recv, task_subtype_bpart_rho);
-          scheduler_activate_recv(s, cj->mpi.recv, task_subtype_bpart_feedback);
+            /* Receive the foreign parts to compute BH accretion rates and do
+             * the swallowing */
+            scheduler_activate_recv(s, cj->mpi.recv, task_subtype_rho);
+            scheduler_activate_recv(s, cj->mpi.recv, task_subtype_part_swallow);
+            scheduler_activate_recv(s, cj->mpi.recv, task_subtype_bpart_merger);
 
-          /* Send the local part information */
-          scheduler_activate_send(s, ci->mpi.send, task_subtype_rho, cj_nodeID);
-          scheduler_activate_send(s, ci->mpi.send, task_subtype_part_swallow,
-                                  cj_nodeID);
-          scheduler_activate_send(s, ci->mpi.send, task_subtype_bpart_merger,
-                                  cj_nodeID);
+            /* Send the local BHs to do feedback */
+            scheduler_activate_send(s, ci->mpi.send,
+                                    task_subtype_bpart_feedback, cj_nodeID);
 
-          /* Drift the cell which will be sent; note that not all sent
-             particles will be drifted, only those that are needed. */
-          cell_activate_drift_part(ci, s);
+            /* Drift before you send */
+            cell_activate_drift_bpart(ci, s);
+          }
+
+          if (cj_active_black_holes) {
+
+            /* Receive the foreign BHs for feedback */
+            scheduler_activate_recv(s, cj->mpi.recv,
+                                    task_subtype_bpart_feedback);
+
+            /* Send the local part information */
+            scheduler_activate_send(s, ci->mpi.send, task_subtype_rho,
+                                    cj_nodeID);
+            scheduler_activate_send(s, ci->mpi.send, task_subtype_part_swallow,
+                                    cj_nodeID);
+            scheduler_activate_send(s, ci->mpi.send, task_subtype_bpart_merger,
+                                    cj_nodeID);
+
+            /* Drift the cell which will be sent; note that not all sent
+               particles will be drifted, only those that are needed. */
+            if (ci->hydro.count > 0) cell_activate_drift_part(ci, s);
+          }
         }
 #endif
       }
@@ -1325,13 +1375,13 @@ void engine_marktasks_mapper(void *map_data, int num_elements,
           /* Is the foreign cell active and will need stuff from us? */
           if (ci_active_gravity) {
 
-            struct link *l = scheduler_activate_send(
-                s, cj->mpi.send, task_subtype_gpart, ci_nodeID);
+            scheduler_activate_send(s, cj->mpi.send, task_subtype_gpart,
+                                    ci_nodeID);
 
             /* Drift the cell which will be sent at the level at which it is
                sent, i.e. drift the cell specified in the send task (l->t)
                itself. */
-            cell_activate_drift_gpart(l->t->ci, s);
+            cell_activate_drift_gpart(cj, s);
           }
 
         } else if (cj_nodeID != nodeID) {
@@ -1343,18 +1393,19 @@ void engine_marktasks_mapper(void *map_data, int num_elements,
           /* Is the foreign cell active and will need stuff from us? */
           if (cj_active_gravity) {
 
-            struct link *l = scheduler_activate_send(
-                s, ci->mpi.send, task_subtype_gpart, cj_nodeID);
+            scheduler_activate_send(s, ci->mpi.send, task_subtype_gpart,
+                                    cj_nodeID);
 
             /* Drift the cell which will be sent at the level at which it is
                sent, i.e. drift the cell specified in the send task (l->t)
                itself. */
-            cell_activate_drift_gpart(l->t->ci, s);
+            cell_activate_drift_gpart(ci, s);
           }
         }
 #endif
-      } /* Only interested in RT tasks as of here. */
+      }
 
+      /* Only interested in RT tasks as of here. */
       else if (t->subtype == task_subtype_rt_gradient) {
 
 #ifdef WITH_MPI
@@ -1527,6 +1578,14 @@ void engine_marktasks_mapper(void *map_data, int num_elements,
         scheduler_activate(s, t);
     }
 
+    /* Star drift tasks? */
+    else if (t_type == task_type_drift_spart) {
+      if (cell_need_activating_stars(t->ci, e, with_star_formation,
+                                     with_star_formation_sink))
+        scheduler_activate(s, t);
+
+    }
+
     /* Star ghost tasks ? */
     else if (t_type == task_type_stars_ghost ||
              t_type == task_type_stars_prep_ghost1 ||
@@ -1542,6 +1601,7 @@ void engine_marktasks_mapper(void *map_data, int num_elements,
       if (cell_need_activating_stars(t->ci, e, with_star_formation,
                                      with_star_formation_sink))
         scheduler_activate(s, t);
+
     }
 
     /* Sink implicit tasks? */
@@ -1556,12 +1616,17 @@ void engine_marktasks_mapper(void *map_data, int num_elements,
              t_type == task_type_bh_swallow_ghost1 ||
              t_type == task_type_bh_swallow_ghost2 ||
              t_type == task_type_bh_swallow_ghost3) {
-      if (cell_is_active_black_holes(t->ci, e)) scheduler_activate(s, t);
+      if (cell_is_active_black_holes(t->ci, e)) {
+
+        scheduler_activate(s, t);
+      }
     }
 
     /* Black holes implicit tasks? */
     else if (t_type == task_type_bh_in || t_type == task_type_bh_out) {
-      if (cell_is_active_black_holes(t->ci, e)) scheduler_activate(s, t);
+      if (cell_is_active_black_holes(t->ci, e)) {
+        scheduler_activate(s, t);
+      }
     }
 
     /* Time-step collection? */
@@ -1631,10 +1696,7 @@ void engine_marktasks_mapper(void *map_data, int num_elements,
 
     /* Radiative transfer implicit tasks */
     else if (t->type == task_type_rt_in) {
-      if (cell_is_rt_active(t->ci, e) ||
-          cell_need_activating_stars(t->ci, e, with_star_formation,
-                                     with_star_formation_sink))
-        scheduler_activate(s, t);
+      if (cell_is_rt_active(t->ci, e)) scheduler_activate(s, t);
     }
 
     else if (t->type == task_type_rt_ghost1 || t->type == task_type_rt_ghost2 ||
@@ -1665,6 +1727,8 @@ void engine_marktasks_mapper(void *map_data, int num_elements,
  */
 int engine_marktasks(struct engine *e) {
 
+#ifdef SWIFT_DEBUG_CHECKS
+
   struct scheduler *s = &e->sched;
   const ticks tic = getticks();
   int rebuild_space = 0;
@@ -1681,4 +1745,9 @@ int engine_marktasks(struct engine *e) {
 
   /* All is well... */
   return rebuild_space;
+
+#else
+  error("Marktasks has been deprecated. Only exists for debugging checks now!");
+  return 0;
+#endif
 }
diff --git a/src/engine_unskip.c b/src/engine_unskip.c
index 30fc66aa16a7e87911a16e149aa6d1c64c1793ce..43af8b5aed2d612f355e29c4134d06e1529fe227 100644
--- a/src/engine_unskip.c
+++ b/src/engine_unskip.c
@@ -160,12 +160,6 @@ static void engine_do_unskip_black_holes(struct cell *c, struct engine *e) {
   /* Early abort (are we below the level where tasks are)? */
   if (!cell_get_flag(c, cell_flag_has_tasks)) return;
 
-  /* Ignore empty cells. */
-  if (c->black_holes.count == 0) return;
-
-  /* Skip inactive cells. */
-  if (!cell_is_active_black_holes(c, e)) return;
-
   /* Recurse */
   if (c->split) {
     for (int k = 0; k < 8; k++) {
diff --git a/src/feedback/EAGLE_thermal/feedback.h b/src/feedback/EAGLE_thermal/feedback.h
index 27bb0c78f13d9d256537f0bc0e081d17736f8e08..aed0ab73934e46577e05953264ae075512f64dee 100644
--- a/src/feedback/EAGLE_thermal/feedback.h
+++ b/src/feedback/EAGLE_thermal/feedback.h
@@ -208,7 +208,8 @@ __attribute__((always_inline)) INLINE static void feedback_prepare_feedback(
     const int with_cosmology) {
 
 #ifdef SWIFT_DEBUG_CHECKS
-  if (sp->birth_time == -1.) error("Evolving a star particle that should not!");
+  if (sp->birth_time == -1. && dt > 0.)
+    error("Evolving a star particle that should not!");
 #endif
 
   /* Start by finishing the loops over neighbours */
diff --git a/src/scheduler.c b/src/scheduler.c
index e0cfee8e102ba731c435c2e8520a93ef9e0d30da..6a9fb140c1f24a66cd029330dc94cffc67ea40b0 100644
--- a/src/scheduler.c
+++ b/src/scheduler.c
@@ -57,6 +57,10 @@
 #include "timers.h"
 #include "version.h"
 
+#ifdef SWIFT_DEBUG_CHECKS
+int activate_by_unskip = 1;
+#endif
+
 /**
  * @brief Re-set the list of active tasks.
  */
@@ -1733,6 +1737,10 @@ struct task *scheduler_addtask(struct scheduler *s, enum task_types type,
   t->tic = 0;
   t->toc = 0;
   t->total_ticks = 0;
+#ifdef SWIFT_DEBUG_CHECKS
+  t->activated_by_unskip = 0;
+  t->activated_by_marktask = 0;
+#endif
 
   if (ci != NULL) cell_set_flag(ci, cell_flag_has_tasks);
   if (cj != NULL) cell_set_flag(cj, cell_flag_has_tasks);
diff --git a/src/scheduler.h b/src/scheduler.h
index b8ec1bfd98ed716d6846f2feec1905b2793f2bab..eb332bf8d1ff5fa655d7af20dd5be1b5ae6fbcdf 100644
--- a/src/scheduler.h
+++ b/src/scheduler.h
@@ -54,6 +54,10 @@
 #define scheduler_flag_none 0
 #define scheduler_flag_steal (1 << 1)
 
+#ifdef SWIFT_DEBUG_CHECKS
+extern int activate_by_unskip;
+#endif
+
 /* Data of a scheduler. */
 struct scheduler {
   /* Scheduler flags. */
@@ -155,6 +159,12 @@ __attribute__((always_inline)) INLINE static void scheduler_activate(
     int ind = atomic_inc(&s->active_count);
     s->tid_active[ind] = t - s->tasks;
   }
+#ifdef SWIFT_DEBUG_CHECKS
+  if (activate_by_unskip)
+    t->activated_by_unskip = 1;
+  else
+    t->activated_by_marktask = 1;
+#endif
 }
 
 /**
diff --git a/src/space_recycle.c b/src/space_recycle.c
index 657d3e23bb1743e8bacc79031bc175174fbf4eb1..cf842273027b22e0123f6efa7a937797e9797f26 100644
--- a/src/space_recycle.c
+++ b/src/space_recycle.c
@@ -142,9 +142,9 @@ void space_rebuild_recycle_mapper(void *map_data, int num_elements,
     c->sinks.do_sink_swallow = NULL;
     c->sinks.do_gas_swallow = NULL;
     c->black_holes.density_ghost = NULL;
-    c->black_holes.swallow_ghost_0 = NULL;
     c->black_holes.swallow_ghost_1 = NULL;
     c->black_holes.swallow_ghost_2 = NULL;
+    c->black_holes.swallow_ghost_3 = NULL;
     c->black_holes.density = NULL;
     c->black_holes.swallow = NULL;
     c->black_holes.do_gas_swallow = NULL;
diff --git a/src/task.h b/src/task.h
index 47649d621b23f0788bf35cd540faaadc287bfe86..b405a0795fd4e00bbc22bf226fe2df9be48335d8 100644
--- a/src/task.h
+++ b/src/task.h
@@ -292,6 +292,12 @@ struct task {
 #ifdef SWIFT_DEBUG_CHECKS
   /* When was this task last run? */
   integertime_t ti_run;
+
+  /* Was this task activted by unskip? */
+  char activated_by_unskip;
+
+  /* Was this task activted by marktask? */
+  char activated_by_marktask;
 #endif /* SWIFT_DEBUG_CHECKS */
 
 } SWIFT_STRUCT_ALIGN;