shared_mem_handler.cc 22.4 KB
Newer Older
Volker Springel's avatar
Volker Springel committed
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
/*******************************************************************************
 * \copyright   This file is part of the GADGET4 N-body/SPH code developed
 * \copyright   by Volker Springel. Copyright (C) 2014-2020 by Volker Springel
 * \copyright   (vspringel@mpa-garching.mpg.de) and all contributing authors.
 *******************************************************************************/

/*! \file  shared_mem_handler.cc
 *
 *  \brief implements code for the shared-memory fetching of remote date through  designated MPI handler ranks
 */

#include <hdf5.h>
#include <mpi.h>
#include <stdio.h>
#include <stdlib.h>
#include <algorithm>
#include <cstring>

#include "../gravtree/gravtree.h"
#include "../ngbtree/ngbtree.h"
#include "../time_integration/driftfac.h"

#include "../mpi_utils/shared_mem_handler.h"

typedef gravtree<simparticles> gtree;
typedef ngbtree ntree;

void shmem::prepare_offset_table(void *p, ptrdiff_t *&offset_tab)  // called by ghost task
{
  ptrdiff_t off = ((char *)p - Mem.Base);

  offset_tab = (ptrdiff_t *)Mem.mymalloc("offset_tab", Island_NTask * sizeof(ptrdiff_t));

  MPI_Gather(&off, sizeof(ptrdiff_t), MPI_BYTE, offset_tab, sizeof(ptrdiff_t), MPI_BYTE, Island_NTask - 1, SharedMemComm);
}

void shmem::inform_offset_table(void *p)  // called by worked tasks
{
  ptrdiff_t off = ((char *)p - Mem.Base);

  MPI_Gather(&off, sizeof(ptrdiff_t), MPI_BYTE, NULL, sizeof(ptrdiff_t), MPI_BYTE, Island_NTask - 1, SharedMemComm);
}

void shmem::free_offset_table(ptrdiff_t *&offset_tab) { Mem.myfree(offset_tab); }

void shmem::shared_memory_handler(void)
{
  simparticles Dp{MPI_COMM_WORLD}; /* dummy needed to access drift functions for ngbtree */

  /* first, we wait for the parameter All.MaxMemSize, so that we can initialize the memory handler */
  MPI_Bcast(All.get_data_ptr(), All.get_data_size(), MPI_BYTE, 0, MPI_COMM_WORLD);
  Mem.mymalloc_init(All.MaxMemSize, RST_BEGIN);

  while(true)
    {
      /* wait for an incoming message */
      MPI_Status status;
      MPI_Probe(MPI_ANY_SOURCE, MPI_ANY_TAG, MPI_COMM_WORLD, &status);

      int source = status.MPI_SOURCE;
      int tag    = status.MPI_TAG;

      int length;
      MPI_Get_count(&status, MPI_BYTE, &length);

      /* now pick it up */
      char *message = (char *)Mem.mymalloc("message", length);
      MPI_Recv(message, length, MPI_BYTE, source, tag, MPI_COMM_WORLD, MPI_STATUS_IGNORE);

      if(tag == TAG_METDATA)  // signals that we are synchronizing addresses and values for tree access
        {
          int handle = *((int *)message);
          Mem.myfree(message);

          MPI_Recv(&All.Ti_Current, sizeof(All.Ti_Current), MPI_BYTE, source, TAG_METDATA + 1, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
          MPI_Recv(&tree_info[handle].Bd, sizeof(bookkeeping_data), MPI_BYTE, source, TAG_METDATA + 2, MPI_COMM_WORLD,
                   MPI_STATUS_IGNORE);

          intposconvert *convfac = &Dp;
          MPI_Recv(convfac, sizeof(intposconvert), MPI_BYTE, source, TAG_METDATA + 3, MPI_COMM_WORLD, MPI_STATUS_IGNORE);

          prepare_offset_table(NULL, tree_info[handle].TopNodes_offsets);
          prepare_offset_table(NULL, tree_info[handle].Nodes_offsets);
          prepare_offset_table(NULL, tree_info[handle].Nextnode_offsets);
          prepare_offset_table(NULL, tree_info[handle].Points_offsets);
          prepare_offset_table(NULL, tree_info[handle].P_offsets);
          prepare_offset_table(NULL, tree_info[handle].SphP_offsets);
          prepare_offset_table(NULL, tree_info[handle].Foreign_Nodes_offsets);
          prepare_offset_table(NULL, tree_info[handle].Foreign_Points_offsets);
        }
      else if(tag == TAG_HEADER)  // signals that we are freeing addresses we stored for tree access
        {
          int handle = *((int *)message);
          Mem.myfree(message);

          free_offset_table(tree_info[handle].Foreign_Points_offsets);
          free_offset_table(tree_info[handle].Foreign_Nodes_offsets);
          free_offset_table(tree_info[handle].SphP_offsets);
          free_offset_table(tree_info[handle].P_offsets);
          free_offset_table(tree_info[handle].Points_offsets);
          free_offset_table(tree_info[handle].Nextnode_offsets);
          free_offset_table(tree_info[handle].Nodes_offsets);
          free_offset_table(tree_info[handle].TopNodes_offsets);
        }
      else if(tag >= TAG_FETCH_SPH_DENSITY && tag < TAG_FETCH_SPH_DENSITY + MAX_TREE_INFOS)  // fetch from SPH density tree
        {
          int handle = tag - TAG_FETCH_SPH_DENSITY;
          deal_with_sph_node_request(message, length, source, handle, &Dp);
          Mem.myfree(message);
        }
      else if(tag >= TAG_FETCH_SPH_HYDRO && tag < TAG_FETCH_SPH_HYDRO + MAX_TREE_INFOS)  // fetch from SPH hydro tree
        {
          int handle = tag - TAG_FETCH_SPH_HYDRO;
          deal_with_sph_node_request(message, length, source, handle, &Dp);
          Mem.myfree(message);
        }
      else if(tag >= TAG_FETCH_SPH_TREETIMESTEP && tag < TAG_FETCH_SPH_TREETIMESTEP + MAX_TREE_INFOS)  // fetch from SPH timesteptree
        {
          int handle = tag - TAG_FETCH_SPH_TREETIMESTEP;
          deal_with_sph_node_request(message, length, source, handle, &Dp);
          Mem.myfree(message);
        }
      else if(tag >= TAG_FETCH_GRAVTREE && tag < TAG_FETCH_GRAVTREE + MAX_TREE_INFOS)  // fetch from gravity tree
        {
          int handle = tag - TAG_FETCH_GRAVTREE;
          deal_with_gravity_node_request(message, length, source, handle);
          Mem.myfree(message);
        }
      else if(tag == TAG_KEY)  // request to terminate gracefully
        {
          H5Eset_auto(H5E_DEFAULT, NULL, NULL);
          MPI_Finalize();
          exit(0);
        }
      else if(tag == TAG_TABLE_ALLOC)  // take over storage for TableData
        {
          size_t tab_len = *((size_t *)message);
          Mem.myfree(message);

          TableData = (char *)Mem.mymalloc("table", tab_len);
          MPI_Recv(TableData, tab_len, MPI_BYTE, source, TAG_DMOM, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
          ptrdiff_t off = ((char *)TableData - Mem.Base);
          MPI_Bcast(&off, sizeof(ptrdiff_t), MPI_BYTE, Island_ThisTask, SharedMemComm);
        }
      else if(tag == TAG_TABLE_FREE)
        {
          Mem.myfree(message);
          Mem.myfree(TableData);
        }
      else if(tag == TAG_EWALD_ALLOC)  // take over storage for EwaldData
        {
          size_t tab_len = *((size_t *)message);
          Mem.myfree(message);

          EwaldData = (char *)Mem.mymalloc("table", tab_len);
          MPI_Recv(EwaldData, tab_len, MPI_BYTE, source, TAG_DMOM, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
          ptrdiff_t off = ((char *)EwaldData - Mem.Base);
          MPI_Bcast(&off, sizeof(ptrdiff_t), MPI_BYTE, Island_ThisTask, SharedMemComm);
        }
      else if(tag == TAG_TOPNODE_ALLOC)  // take over the storage of common top-level tree data
        {
          int handle = num_tree_info++;
          if(num_tree_info > MAX_TREE_INFOS)
            Terminate("num_tree_info > MAX_TREE_INFOS");

          size_t *sizep   = ((size_t *)message);
          size_t sizes[4] = {sizep[0], sizep[1], sizep[2], sizep[3]};
          Mem.myfree(message);

          tree_info[handle].NodeLevel_storage   = (char *)Mem.mymalloc("NodeLevel_storage", sizes[0]);
          tree_info[handle].NodeSibling_storage = (char *)Mem.mymalloc("NodeSibling_storage", sizes[1]);
          tree_info[handle].NodeIndex_storage   = (char *)Mem.mymalloc("NodeIndex_storage", sizes[2]);
          tree_info[handle].TopNodes_storage    = (char *)Mem.mymalloc("TopNodes_storage", sizes[3]);

          ptrdiff_t off[4] = {
              ((char *)tree_info[handle].NodeLevel_storage - Mem.Base), ((char *)tree_info[handle].NodeSibling_storage - Mem.Base),
              ((char *)tree_info[handle].NodeIndex_storage - Mem.Base), ((char *)tree_info[handle].TopNodes_storage - Mem.Base)};

          MPI_Send(off, 4 * sizeof(ptrdiff_t), MPI_BYTE, source, TAG_TOPNODE_OFFSET, MPI_COMM_WORLD);

          MPI_Send(&handle, 1, MPI_INT, source, TAG_N, MPI_COMM_WORLD);
        }
      else if(tag == TAG_TOPNODE_FREE)  // free the top-level storage for a tree again
        {
          int handle = *((int *)message);
          Mem.myfree(message);

          num_tree_info--;
          if(handle != num_tree_info)
            Terminate("unexpected handle");

          Mem.myfree(tree_info[handle].TopNodes_storage);
          Mem.myfree(tree_info[handle].NodeIndex_storage);
          Mem.myfree(tree_info[handle].NodeSibling_storage);
          Mem.myfree(tree_info[handle].NodeLevel_storage);
        }
      else if(tag == TAG_DRIFT_INIT)  // make the shared memory handler update All and init the local drift tables
        {
          memcpy(All.get_data_ptr(), message, All.get_data_size());
          Mem.myfree(message);
          Driftfac.init_drift_table();
        }
203
204
205
206
207
      else if(tag == TAG_ALL_UPDATE)  // make the shared memory handler update the contents of the All structure
        {
          memcpy(All.get_data_ptr(), message, All.get_data_size());
          Mem.myfree(message);
        }
Volker Springel's avatar
Volker Springel committed
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
    }
}

void shmem::deal_with_sph_node_request(char *message, int length, int source, int handle, simparticles *Sp)
{
#ifndef LEAN
  bookkeeping_data &Bdat = tree_info[handle].Bd;

  // we got the list of requested nodes
  ntree::node_req *node_req_recv = (ntree::node_req *)message;
  int nrecv                      = length / sizeof(ntree::node_req);

  /* as part of this message, we get the real actually targeted rank in
   * the simulation communicator. We will translate this to the rank in our shared memory
   * block
   */

  /* now prepare answer message by reading from shared memory */
  /******************* prepare tree answer ********************/

  ntree::node_count_info *node_info_recv =
      (ntree::node_count_info *)Mem.mymalloc("node_info_recv", nrecv * sizeof(ntree::node_count_info));

  /* first let's count how many nodes and particles are hanging below in each case */
  int n_recvpoints = 0;
  int n_recvnodes  = 0;

  for(int i = 0; i < nrecv; i++)
    {
      node_info_recv[i].count_nodes = 0;
      node_info_recv[i].count_parts = 0;

      int no      = node_req_recv[i].foreignnode;
      int task    = node_req_recv[i].foreigntask;
      int shmrank = GetShmRankForSimulCommRank[task];  // corresponding local shared memory rank, stays fixed below

      if(no < Bdat.MaxPart || no >= Bdat.MaxPart + Bdat.MaxNodes)
        Terminate("not an internal node");

      ngbnode *nop = ((ngbnode *)get_basenodep(no, shmrank, handle)) + no;

      int p = nop->nextnode;

      while(p != nop->sibling)
        {
          if(p < 0)
            Terminate("p=%d < 0", p);

          if(p < Bdat.MaxPart) /* a local particle */
            {
              node_info_recv[i].count_parts++;
              p = get_nextnodep(shmrank, handle)[p];
            }
          else if(p < Bdat.MaxPart + Bdat.MaxNodes) /* an internal node  */
            {
              node_info_recv[i].count_nodes++;
              p = (((ngbnode *)get_basenodep(p, shmrank, handle)) + p)->sibling;
            }
          else if(p >= Bdat.ImportedNodeOffset && p < Bdat.EndOfTreePoints) /* an imported tree point */
            {
              node_info_recv[i].count_parts++;
              p = get_nextnodep(shmrank, handle)[p - Bdat.MaxNodes];
            }
          else
            Terminate("p=%d MaxPart=%d MaxNodes=%d", p, Bdat.MaxPart, Bdat.MaxNodes);
        }

      if(node_info_recv[i].count_parts == 0 && node_info_recv[i].count_nodes == 0)
        Terminate("strange: we have we requested an empty node?\n");

      n_recvpoints += node_info_recv[i].count_parts;
      n_recvnodes += node_info_recv[i].count_nodes;
    }

  foreign_sphpoint_data *exportbuf_points = (foreign_sphpoint_data *)Mem.mymalloc_movable(
      &exportbuf_points, "exportbuf_points", n_recvpoints * sizeof(foreign_sphpoint_data));
  ngbnode *exportbuf_nodes = (ngbnode *)Mem.mymalloc_movable(&exportbuf_nodes, "exportbuf_nodes", n_recvnodes * sizeof(ngbnode));

  n_recvpoints = 0;
  n_recvnodes  = 0;

  for(int i = 0; i < nrecv; i++)
    {
      int no      = node_req_recv[i].foreignnode;
      int task    = node_req_recv[i].foreigntask;
      int shmrank = GetShmRankForSimulCommRank[task];  // corresponding local shared memory rank, stays fixed below

      ngbnode *nop = ((ngbnode *)get_basenodep(no, shmrank, handle)) + no;

      int p = nop->nextnode;

      while(p != nop->sibling)
        {
          if(p < Bdat.MaxPart) /* a local particle */
            {
              int off = n_recvpoints++;

              foreign_sphpoint_data *expoints = &exportbuf_points[off];

              particle_data *ptr         = (particle_data *)get_Pp(shmrank, handle) + p;
              sph_particle_data *sph_ptr = (sph_particle_data *)get_SphPp(shmrank, handle) + p;

              particle_data p_copy;
              sph_particle_data sphp_copy;

              if(ptr->get_Ti_Current() != All.Ti_Current)
                {
                  /* Because of possible lightcone output, shared memory fetch not allowed to drift original,
                   * because this rank doesn't have a lightcone buffer. The original node needs to do it.
                   * We thus drift a copy of the particle without allowing lightcone access
                   */
#ifndef LEAN
                  while(ptr->access.test_and_set(std::memory_order_acquire))
                    ;  // acquire spin lock
#endif
                  p_copy    = *ptr;
                  sphp_copy = *sph_ptr;

#ifndef LEAN
                  ptr->access.clear(std::memory_order_release);  // release spin lock
#endif
                  p_copy.access.clear(std::memory_order_release);  // clear spin lock in copy

                  /* use the copy from now on */

                  ptr     = &p_copy;
                  sph_ptr = &sphp_copy;

                  // the final flag tells the drift to not consider lightcone crossings
                  Sp->drift_particle(ptr, sph_ptr, All.Ti_Current, true);
                }

              expoints->IntPos[0]    = ptr->IntPos[0];
              expoints->IntPos[1]    = ptr->IntPos[1];
              expoints->IntPos[2]    = ptr->IntPos[2];
              expoints->Mass         = ptr->getMass();
              expoints->TimeBinHydro = ptr->TimeBinHydro;
              expoints->SphCore      = *sph_ptr;

              expoints->Nextnode = -1;

              p = get_nextnodep(shmrank, handle)[p];
            }
          else if(p < Bdat.MaxPart + Bdat.MaxNodes) /* an internal node  */
            {
              int off = n_recvnodes++;

              ngbnode *sourcep = (((ngbnode *)get_basenodep(p, shmrank, handle)) + p);
              ngbnode *exnodes = &exportbuf_nodes[off];

              if(sourcep->Ti_Current != All.Ti_Current)
                sourcep->drift_node(All.Ti_Current, Sp);

              *exnodes = *sourcep;

              exnodes->cannot_be_opened_locally = 1;
              exnodes->nextnode                 = -1;
              exnodes->sibling                  = -1;
              exnodes->OriginTask               = task;
              exnodes->OriginNode               = p;

              p = sourcep->sibling;
            }
          else if(p >= Bdat.ImportedNodeOffset) /* an imported treepoint particle  */
            {
              Terminate("not expected here");
            }
        }
    }

  /************************************************************/

  MPI_Send(node_info_recv, nrecv * sizeof(ntree::node_count_info), MPI_BYTE, source, TAG_N, MPI_COMM_WORLD);

  /* now transfer the points and nodes */
  if(n_recvpoints > 0)
    MPI_Send(exportbuf_points, n_recvpoints * sizeof(foreign_sphpoint_data), MPI_BYTE, source, TAG_PDATA, MPI_COMM_WORLD);

  if(n_recvnodes > 0)
    MPI_Send(exportbuf_nodes, n_recvnodes * sizeof(ngbnode), MPI_BYTE, source, TAG_SPHDATA, MPI_COMM_WORLD);

  Mem.myfree(exportbuf_nodes);
  Mem.myfree(exportbuf_points);
  Mem.myfree(node_info_recv);
#endif
}

void shmem::deal_with_gravity_node_request(char *message, int length, int source, int handle)
{
  bookkeeping_data &Bdat = tree_info[handle].Bd;

  // we got the list of requested nodes
  gtree::node_req *node_req_recv = (gtree::node_req *)message;
  int nrecv                      = length / sizeof(gtree::node_req);

  /* as part of this message, we get the real actually targeted rank in
   * the simulation communicator. We will translate this to the rank in our shared memory
   * block
   */

  /* now prepare answer message by reading from shared memory */
  /******************* prepare tree answer ********************/

  gtree::node_count_info *node_info_recv =
      (gtree::node_count_info *)Mem.mymalloc("node_info_recv", nrecv * sizeof(gtree::node_count_info));

  /* first let's count how many nodes and particles are hanging below in each case */
  int n_recvpoints = 0;
  int n_recvnodes  = 0;

  for(int i = 0; i < nrecv; i++)
    {
      node_info_recv[i].count_nodes = 0;
      node_info_recv[i].count_parts = 0;

      int no      = node_req_recv[i].foreignnode;
      int task    = node_req_recv[i].foreigntask;
      int shmrank = GetShmRankForSimulCommRank[task];  // corresponding local shared memory rank, stays fixed below

      if(no < Bdat.MaxPart || no >= Bdat.MaxPart + Bdat.MaxNodes)
        Terminate("not an internal node");

      gravnode *nop = ((gravnode *)get_basenodep(no, shmrank, handle)) + no;

      int p = nop->nextnode;

      while(p != nop->sibling)
        {
          if(p < 0)
            Terminate("p=%d < 0", p);

          if(p < Bdat.MaxPart) /* a local particle */
            {
              node_info_recv[i].count_parts++;
              p = get_nextnodep(shmrank, handle)[p];
            }
          else if(p < Bdat.MaxPart + Bdat.MaxNodes) /* an internal node  */
            {
              node_info_recv[i].count_nodes++;
              p = (((gravnode *)get_basenodep(p, shmrank, handle)) + p)->sibling;
            }
          else if(p >= Bdat.ImportedNodeOffset && p < Bdat.EndOfTreePoints) /* an imported tree point */
            {
              node_info_recv[i].count_parts++;
              p = get_nextnodep(shmrank, handle)[p - Bdat.MaxNodes];
            }
          else
            Terminate("p=%d MaxPart=%d MaxNodes=%d", p, Bdat.MaxPart, Bdat.MaxNodes);
        }

      if(node_info_recv[i].count_parts == 0 && node_info_recv[i].count_nodes == 0)
        Terminate("strange: we have we requested an empty node?\n");

      n_recvpoints += node_info_recv[i].count_parts;
      n_recvnodes += node_info_recv[i].count_nodes;
    }

  foreign_gravpoint_data *exportbuf_points = (foreign_gravpoint_data *)Mem.mymalloc_movable(
      &exportbuf_points, "exportbuf_points", n_recvpoints * sizeof(foreign_gravpoint_data));
  gravnode *exportbuf_nodes = (gravnode *)Mem.mymalloc_movable(&exportbuf_nodes, "exportbuf_nodes", n_recvnodes * sizeof(gravnode));

  n_recvpoints = 0;
  n_recvnodes  = 0;

  for(int i = 0; i < nrecv; i++)
    {
      int no      = node_req_recv[i].foreignnode;
      int task    = node_req_recv[i].foreigntask;
      int shmrank = GetShmRankForSimulCommRank[task];  // corresponding local shared memory rank, stays fixed below

      gravnode *nop = ((gravnode *)get_basenodep(no, shmrank, handle)) + no;

      int p = nop->nextnode;

      while(p != nop->sibling)
        {
          if(p < Bdat.MaxPart) /* a local particle */
            {
              int off = n_recvpoints++;

              foreign_gravpoint_data *expoints = &exportbuf_points[off];

              particle_data *ptr = (particle_data *)get_Pp(shmrank, handle) + p;

              expoints->IntPos[0] = ptr->IntPos[0];
              expoints->IntPos[1] = ptr->IntPos[1];
              expoints->IntPos[2] = ptr->IntPos[2];
              expoints->Mass      = ptr->getMass();
              expoints->Type      = ptr->getType();
              expoints->OldAcc    = ptr->OldAcc;
#if NSOFTCLASSES > 1
              expoints->SofteningClass = ptr->getSofteningClass();
#endif
#if defined(PMGRID) && defined(PLACEHIGHRESREGION)
              expoints->InsideOutsideFlag = ptr->InsideOutsideFlag;
#endif
              expoints->Nextnode = -1;

              p = get_nextnodep(shmrank, handle)[p];
            }
          else if(p < Bdat.MaxPart + Bdat.MaxNodes) /* an internal node  */
            {
              int off = n_recvnodes++;

              gravnode *exnodes = &exportbuf_nodes[off];
              gravnode *sourcep = (((gravnode *)get_basenodep(p, shmrank, handle)) + p);

              memcpy(static_cast<void *>(exnodes), static_cast<void *>(sourcep),
                     sizeof(gravnode));  //  cannot do a  *exnodes = *sourcep; because out std::atomic_flag
                                         //  has a deleted default copy operator

              exnodes->cannot_be_opened_locally = 1;
              exnodes->nextnode                 = -1;
              exnodes->sibling                  = -1;
              exnodes->OriginTask               = task;
              exnodes->OriginNode               = p;

              p = sourcep->sibling;
            }
          else if(p >= Bdat.ImportedNodeOffset) /* an imported Treepoint particle  */
            {
              int off = n_recvpoints++;

              foreign_gravpoint_data *expoints = &exportbuf_points[off];
              int n                            = p - Bdat.ImportedNodeOffset;
              gravpoint_data *pointsp          = ((gravpoint_data *)get_pointsp(shmrank, handle)) + n;

              expoints->IntPos[0] = pointsp->IntPos[0];
              expoints->IntPos[1] = pointsp->IntPos[1];
              expoints->IntPos[2] = pointsp->IntPos[2];
              expoints->Mass      = pointsp->Mass;
              expoints->Type      = pointsp->Type;
              expoints->OldAcc    = pointsp->OldAcc;
#if NSOFTCLASSES > 1
              expoints->SofteningClass = pointsp->SofteningClass;
#endif
#if defined(PMGRID) && defined(PLACEHIGHRESREGION)
              expoints->InsideOutsideFlag = pointsp->InsideOutsideFlag;
#endif
              expoints->Nextnode = -1;

              p = get_nextnodep(shmrank, handle)[p - Bdat.MaxNodes];
            }
        }
    }

  /************************************************************/

  MPI_Send(node_info_recv, nrecv * sizeof(gtree::node_count_info), MPI_BYTE, source, TAG_N, MPI_COMM_WORLD);

  /* now transfer the points and nodes */
  if(n_recvpoints > 0)
    MPI_Send(exportbuf_points, n_recvpoints * sizeof(foreign_gravpoint_data), MPI_BYTE, source, TAG_PDATA, MPI_COMM_WORLD);

  if(n_recvnodes > 0)
    MPI_Send(exportbuf_nodes, n_recvnodes * sizeof(gravnode), MPI_BYTE, source, TAG_SPHDATA, MPI_COMM_WORLD);

  Mem.myfree(exportbuf_nodes);
  Mem.myfree(exportbuf_points);
  Mem.myfree(node_info_recv);
}