19#include "../gravtree/gravtree.h"
20#include "../ngbtree/ngbtree.h"
21#include "../time_integration/driftfac.h"
23#include "../mpi_utils/shared_mem_handler.h"
30 ptrdiff_t off = ((
char *)p -
Mem.
Base);
32 offset_tab = (ptrdiff_t *)
Mem.mymalloc(
"offset_tab",
Island_NTask *
sizeof(ptrdiff_t));
39 ptrdiff_t off = ((
char *)p -
Mem.
Base);
58 MPI_Probe(MPI_ANY_SOURCE, MPI_ANY_TAG, MPI_COMM_WORLD, &status);
60 int source = status.MPI_SOURCE;
61 int tag = status.MPI_TAG;
64 MPI_Get_count(&status, MPI_BYTE, &length);
67 char *message = (
char *)
Mem.mymalloc(
"message", length);
68 MPI_Recv(message, length, MPI_BYTE, source, tag, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
72 int handle = *((
int *)message);
96 int handle = *((
int *)message);
134 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
140 size_t tab_len = *((
size_t *)message);
144 MPI_Recv(
TableData, tab_len, MPI_BYTE, source,
TAG_DMOM, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
155 size_t tab_len = *((
size_t *)message);
159 MPI_Recv(
EwaldData, tab_len, MPI_BYTE, source,
TAG_DMOM, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
167 Terminate(
"num_tree_info > MAX_TREE_INFOS");
169 size_t *sizep = ((
size_t *)message);
170 size_t sizes[4] = {sizep[0], sizep[1], sizep[2], sizep[3]};
182 MPI_Send(off, 4 *
sizeof(ptrdiff_t), MPI_BYTE, source,
TAG_TOPNODE_OFFSET, MPI_COMM_WORLD);
184 MPI_Send(&handle, 1, MPI_INT, source,
TAG_N, MPI_COMM_WORLD);
188 int handle = *((
int *)message);
220 ntree::node_req *node_req_recv = (ntree::node_req *)message;
221 int nrecv = length /
sizeof(ntree::node_req);
231 ntree::node_count_info *node_info_recv =
232 (ntree::node_count_info *)
Mem.mymalloc(
"node_info_recv", nrecv *
sizeof(ntree::node_count_info));
235 int n_recvpoints = 0;
238 for(
int i = 0; i < nrecv; i++)
240 node_info_recv[i].count_nodes = 0;
241 node_info_recv[i].count_parts = 0;
243 int no = node_req_recv[i].foreignnode;
244 int task = node_req_recv[i].foreigntask;
261 node_info_recv[i].count_parts++;
266 node_info_recv[i].count_nodes++;
271 node_info_recv[i].count_parts++;
278 if(node_info_recv[i].count_parts == 0 && node_info_recv[i].count_nodes == 0)
279 Terminate(
"strange: we have we requested an empty node?\n");
281 n_recvpoints += node_info_recv[i].count_parts;
282 n_recvnodes += node_info_recv[i].count_nodes;
287 ngbnode *exportbuf_nodes = (
ngbnode *)
Mem.mymalloc_movable(&exportbuf_nodes,
"exportbuf_nodes", n_recvnodes *
sizeof(
ngbnode));
292 for(
int i = 0; i < nrecv; i++)
294 int no = node_req_recv[i].foreignnode;
295 int task = node_req_recv[i].foreigntask;
306 int off = n_recvpoints++;
323 while(ptr->
access.test_and_set(std::memory_order_acquire))
327 sphp_copy = *sph_ptr;
330 ptr->
access.clear(std::memory_order_release);
332 p_copy.
access.clear(std::memory_order_release);
337 sph_ptr = &sphp_copy;
356 int off = n_recvnodes++;
359 ngbnode *exnodes = &exportbuf_nodes[off];
383 MPI_Send(node_info_recv, nrecv *
sizeof(ntree::node_count_info), MPI_BYTE, source,
TAG_N, MPI_COMM_WORLD);
390 MPI_Send(exportbuf_nodes, n_recvnodes *
sizeof(
ngbnode), MPI_BYTE, source,
TAG_SPHDATA, MPI_COMM_WORLD);
392 Mem.myfree(exportbuf_nodes);
393 Mem.myfree(exportbuf_points);
394 Mem.myfree(node_info_recv);
403 gtree::node_req *node_req_recv = (gtree::node_req *)message;
404 int nrecv = length /
sizeof(gtree::node_req);
414 gtree::node_count_info *node_info_recv =
415 (gtree::node_count_info *)
Mem.mymalloc(
"node_info_recv", nrecv *
sizeof(gtree::node_count_info));
418 int n_recvpoints = 0;
421 for(
int i = 0; i < nrecv; i++)
423 node_info_recv[i].count_nodes = 0;
424 node_info_recv[i].count_parts = 0;
426 int no = node_req_recv[i].foreignnode;
427 int task = node_req_recv[i].foreigntask;
444 node_info_recv[i].count_parts++;
449 node_info_recv[i].count_nodes++;
454 node_info_recv[i].count_parts++;
461 if(node_info_recv[i].count_parts == 0 && node_info_recv[i].count_nodes == 0)
462 Terminate(
"strange: we have we requested an empty node?\n");
464 n_recvpoints += node_info_recv[i].count_parts;
465 n_recvnodes += node_info_recv[i].count_nodes;
475 for(
int i = 0; i < nrecv; i++)
477 int no = node_req_recv[i].foreignnode;
478 int task = node_req_recv[i].foreigntask;
489 int off = n_recvpoints++;
504#if defined(PMGRID) && defined(PLACEHIGHRESREGION)
505 expoints->InsideOutsideFlag = ptr->InsideOutsideFlag;
513 int off = n_recvnodes++;
515 gravnode *exnodes = &exportbuf_nodes[off];
518 memcpy(
static_cast<void *
>(exnodes),
static_cast<void *
>(sourcep),
532 int off = n_recvpoints++;
547#if defined(PMGRID) && defined(PLACEHIGHRESREGION)
548 expoints->InsideOutsideFlag = pointsp->InsideOutsideFlag;
559 MPI_Send(node_info_recv, nrecv *
sizeof(gtree::node_count_info), MPI_BYTE, source,
TAG_N, MPI_COMM_WORLD);
566 MPI_Send(exportbuf_nodes, n_recvnodes *
sizeof(
gravnode), MPI_BYTE, source,
TAG_SPHDATA, MPI_COMM_WORLD);
568 Mem.myfree(exportbuf_nodes);
569 Mem.myfree(exportbuf_points);
570 Mem.myfree(node_info_recv);
global_data_all_processes All
void init_drift_table(void)
void mymalloc_init(int maxmemsize, enum restart_options restartflag)
Initialize memory manager.
void prepare_offset_table(void *p, ptrdiff_t *&offset_tab)
void deal_with_gravity_node_request(char *message, int length, int source, int handle)
int * get_nextnodep(unsigned char shmrank, int handle)
void deal_with_sph_node_request(char *message, int length, int source, int handle, simparticles *Sp)
void shared_memory_handler(void)
void inform_offset_table(void *p)
int * GetShmRankForSimulCommRank
tree_storage_info tree_info[MAX_TREE_INFOS]
char * get_basenodep(int no, unsigned char shmrank, int handle)
char * get_pointsp(unsigned char shmrank, int handle)
char * get_Pp(unsigned char shmrank, int handle)
char * get_SphPp(unsigned char shmrank, int handle)
void free_offset_table(ptrdiff_t *&offset_tab)
int drift_particle(particle_data *P, sph_particle_data *SphP, integertime time1, bool ignore_light_cone=false)
This function drifts a particle i to time1.
#define TAG_FETCH_SPH_DENSITY
#define TAG_FETCH_SPH_TREETIMESTEP
#define TAG_FETCH_GRAVTREE
#define TAG_TOPNODE_ALLOC
#define TAG_FETCH_SPH_HYDRO
#define TAG_TOPNODE_OFFSET
gravtree< simparticles > gtree
std::atomic< unsigned char > cannot_be_opened_locally
unsigned char SofteningClass
sph_particle_data_hydrocore SphCore
size_t get_data_size(void)
char * get_data_ptr(void)
unsigned char SofteningClass
std::atomic< integertime > Ti_Current
void drift_node(integertime time1, simparticles *Sp)
integertime get_Ti_Current(void)
unsigned char getSofteningClass(void)
unsigned char getType(void)
char * NodeSibling_storage