Particle *move_indexed_particle(ParticleList *dl, ParticleList *sl, int i) { int re = realloc_particlelist(dl, ++dl->n); Particle *dst = &dl->part[dl->n - 1]; Particle *src = &sl->part[i]; Particle *end = &sl->part[sl->n - 1]; memcpy(dst, src, sizeof(Particle)); if (re) { //fprintf(stderr, "%d: m_i_p: update destination list after realloc\n",this_node); update_local_particles(dl); } else { //fprintf(stderr, "%d: m_i_p: update loc_part entry for moved particle (id %d)\n",this_node,dst->p.identity); local_particles[dst->p.identity] = dst; } if ( src != end ) { //fprintf(stderr, "%d: m_i_p: copy end particle in source list (id %d)\n",this_node,end->p.identity); memcpy(src, end, sizeof(Particle)); } if (realloc_particlelist(sl, --sl->n)) { //fprintf(stderr, "%d: m_i_p: update source list after realloc\n",this_node); update_local_particles(sl); } else if ( src != end ) { //fprintf(stderr, "%d: m_i_p: update loc_part entry for end particle (id %d)\n",this_node,src->p.identity); local_particles[src->p.identity] = src; } return dst; }
void local_sort_particles() { CELL_TRACE(fprintf(stderr, "%d: entering local_sort_particles\n", this_node)); /* first distribute strictly on nodes */ cells_resort_particles(CELL_GLOBAL_EXCHANGE); CELL_TRACE(fprintf(stderr, "%d: sorting local cells\n", this_node)); /* now sort the local cells */ for (int c = 0; c < local_cells.n; c++) { Cell *cell = local_cells.cell[c]; Particle *p = cell->part; int np = cell->n; #ifdef CELL_DEBUG for (int id = 0; id < np; ++id) { Cell *tgt_cell = cell_structure.position_to_cell(p[id].r.p); if (tgt_cell != cell) { fprintf(stderr, "%d: particle %d at position %lf %lf %lf is not in its expected cell. Have %ld, expected %ld\n", this_node, p[id].p.identity, p[id].r.p[0], p[id].r.p[1], p[id].r.p[2], (cell - cells)/sizeof(Cell*), (tgt_cell - cells) /sizeof(Cell*)); } } #endif qsort(p, np, sizeof(Particle), compare_particles); update_local_particles(cell); } CELL_TRACE(dump_particle_ordering()); CELL_TRACE(fprintf(stderr, "%d: leaving local_sort_particles\n", this_node)); }
void local_place_particle(int part, double p[3], int _new) { Cell *cell; double pp[3]; int i[3], rl; Particle *pt; i[0] = 0; i[1] = 0; i[2] = 0; pp[0] = p[0]; pp[1] = p[1]; pp[2] = p[2]; #ifdef LEES_EDWARDS double vv[3]={0.,0.,0.}; fold_position(pp, vv, i); #else fold_position(pp, i); #endif if (_new) { /* allocate particle anew */ cell = cell_structure.position_to_cell(pp); if (!cell) { fprintf(stderr, "%d: INTERNAL ERROR: particle %d at %f(%f) %f(%f) %f(%f) does not belong on this node\n", this_node, part, p[0], pp[0], p[1], pp[1], p[2], pp[2]); errexit(); } rl = realloc_particlelist(cell, ++cell->n); pt = &cell->part[cell->n - 1]; init_particle(pt); pt->p.identity = part; if (rl) update_local_particles(cell); else local_particles[pt->p.identity] = pt; } else pt = local_particles[part]; PART_TRACE(fprintf(stderr, "%d: local_place_particle: got particle id=%d @ %f %f %f\n", this_node, part, p[0], p[1], p[2])); #ifdef LEES_EDWARDS pt->m.v[0] += vv[0]; pt->m.v[1] += vv[1]; pt->m.v[2] += vv[2]; #endif memcpy(pt->r.p, pp, 3*sizeof(double)); memcpy(pt->l.i, i, 3*sizeof(int)); #ifdef BOND_CONSTRAINT memcpy(pt->r.p_old, pp, 3*sizeof(double)); #endif }
Particle *append_indexed_particle(ParticleList *l, Particle *part) { int re; Particle *p; re = realloc_particlelist(l, ++l->n); p = &l->part[l->n - 1]; memcpy(p, part, sizeof(Particle)); if (re) update_local_particles(l); else local_particles[p->p.identity] = p; return p; }
void nsq_topology_init(CellPList *old) { Particle *part; int n, c, p, np, ntodo, diff; CELL_TRACE(fprintf(stderr, "%d: nsq_topology_init, %d\n", this_node, old->n)); cell_structure.type = CELL_STRUCTURE_NSQUARE; cell_structure.position_to_node = map_position_node_array; cell_structure.position_to_cell = nsq_position_to_cell; realloc_cells(n_nodes); /* mark cells */ local = &cells[this_node]; realloc_cellplist(&local_cells, local_cells.n = 1); local_cells.cell[0] = local; realloc_cellplist(&ghost_cells, ghost_cells.n = n_nodes - 1); c = 0; for (n = 0; n < n_nodes; n++) if (n != this_node) ghost_cells.cell[c++] = &cells[n]; /* distribute force calculation work */ ntodo = (n_nodes + 3)/2; init_cellplist(&me_do_ghosts); realloc_cellplist(&me_do_ghosts, ntodo); for (n = 0; n < n_nodes; n++) { diff = n - this_node; /* simple load balancing formula. Basically diff % n, where n >= n_nodes, n odd. The node itself is also left out, as it is treated differently */ if (((diff > 0 && (diff % 2) == 0) || (diff < 0 && ((-diff) % 2) == 1))) { CELL_TRACE(fprintf(stderr, "%d: doing interactions with %d\n", this_node, n)); me_do_ghosts.cell[me_do_ghosts.n++] = &cells[n]; } } /* create communicators */ nsq_prepare_comm(&cell_structure.ghost_cells_comm, GHOSTTRANS_PARTNUM); nsq_prepare_comm(&cell_structure.exchange_ghosts_comm, GHOSTTRANS_PROPRTS | GHOSTTRANS_POSITION); nsq_prepare_comm(&cell_structure.update_ghost_pos_comm, GHOSTTRANS_POSITION); nsq_prepare_comm(&cell_structure.collect_ghost_force_comm, GHOSTTRANS_FORCE); /* here we just decide what to transfer where */ if (n_nodes > 1) { for (n = 0; n < n_nodes; n++) { /* use the prefetched send buffers. Node 0 transmits first and never prefetches. */ if (this_node == 0 || this_node != n) { cell_structure.ghost_cells_comm.comm[n].type = GHOST_BCST; cell_structure.exchange_ghosts_comm.comm[n].type = GHOST_BCST; cell_structure.update_ghost_pos_comm.comm[n].type = GHOST_BCST; } else { cell_structure.ghost_cells_comm.comm[n].type = GHOST_BCST | GHOST_PREFETCH; cell_structure.exchange_ghosts_comm.comm[n].type = GHOST_BCST | GHOST_PREFETCH; cell_structure.update_ghost_pos_comm.comm[n].type = GHOST_BCST | GHOST_PREFETCH; } cell_structure.collect_ghost_force_comm.comm[n].type = GHOST_RDCE; } /* first round: all nodes except the first one prefetch their send data */ if (this_node != 0) { cell_structure.ghost_cells_comm.comm[0].type |= GHOST_PREFETCH; cell_structure.exchange_ghosts_comm.comm[0].type |= GHOST_PREFETCH; cell_structure.update_ghost_pos_comm.comm[0].type |= GHOST_PREFETCH; } } /* copy particles */ for (c = 0; c < old->n; c++) { part = old->cell[c]->part; np = old->cell[c]->n; for (p = 0; p < np; p++) append_unindexed_particle(local, &part[p]); } update_local_particles(local); }
void nsq_balance_particles() { int i, n, surplus, s_node, tmp, lack, l_node, transfer; int pp = cells_get_n_particles(); int *ppnode = malloc(n_nodes*sizeof(int)); /* minimal difference between node shares */ int minshare = n_total_particles/n_nodes; int maxshare = minshare + 1; CELL_TRACE(fprintf(stderr, "%d: nsq_balance_particles: load %d-%d\n", this_node, minshare, maxshare)); MPI_Allgather(&pp, 1, MPI_INT, ppnode, 1, MPI_INT, MPI_COMM_WORLD); for (;;) { /* find node with most excessive particles */ surplus = -1; s_node = -1; for (n = 0; n < n_nodes; n++) { tmp = ppnode[n] - minshare; CELL_TRACE(fprintf(stderr, "%d: nsq_balance_particles: node %d has %d\n", this_node, n, ppnode[n])); if (tmp > surplus) { surplus = tmp; s_node = n; } } CELL_TRACE(fprintf(stderr, "%d: nsq_balance_particles: excess %d on node %d\n", this_node, surplus, s_node)); /* find node with most lacking particles */ lack = -1; l_node = -1; for (n = 0; n < n_nodes; n++) { tmp = maxshare - ppnode[n]; if (tmp > lack) { lack = tmp; l_node = n; } } CELL_TRACE(fprintf(stderr, "%d: nsq_balance_particles: lack %d on node %d\n", this_node, lack, l_node)); /* should not happen: minshare or maxshare wrong or more likely, the algorithm */ if (s_node == -1 || l_node == -1) { fprintf(stderr, "%d: Particle load balancing failed\n", this_node); break; } /* exit if all nodes load is withing min and max share */ if (lack <= 1 && surplus <= 1) break; transfer = lack < surplus ? lack : surplus; if (s_node == this_node) { ParticleList send_buf; init_particlelist(&send_buf); realloc_particlelist(&send_buf, send_buf.n = transfer); for (i = 0; i < transfer; i++) { memcpy(&send_buf.part[i], &local->part[--local->n], sizeof(Particle)); } realloc_particlelist(local, local->n); update_local_particles(local); send_particles(&send_buf, l_node); #ifdef ADDITIONAL_CHECKS check_particle_consistency(); #endif } else if (l_node == this_node) { recv_particles(local, s_node); #ifdef ADDITIONAL_CHECKS check_particle_consistency(); #endif } ppnode[s_node] -= transfer; ppnode[l_node] += transfer; } CELL_TRACE(fprintf(stderr, "%d: nsq_balance_particles: done\n", this_node)); free(ppnode); }
void recv_particles(ParticleList *particles, int node) { int transfer=0, read, pc; IntList local_dyn; PART_TRACE(fprintf(stderr, "%d: recv_particles from %d\n", this_node, node)); MPI_Recv(&transfer, 1, MPI_INT, node, REQ_SNDRCV_PART, comm_cart, MPI_STATUS_IGNORE); PART_TRACE(fprintf(stderr, "%d: recv_particles get %d\n", this_node, transfer)); realloc_particlelist(particles, particles->n + transfer); MPI_Recv(&particles->part[particles->n], transfer*sizeof(Particle), MPI_BYTE, node, REQ_SNDRCV_PART, comm_cart, MPI_STATUS_IGNORE); particles->n += transfer; init_intlist(&local_dyn); for (pc = particles->n - transfer; pc < particles->n; pc++) { Particle *p = &particles->part[pc]; local_dyn.n += p->bl.n; #ifdef EXCLUSIONS local_dyn.n += p->el.n; #endif PART_TRACE(fprintf(stderr, "%d: recv_particles got particle %d\n", this_node, p->p.identity)); if (local_particles[p->p.identity] != NULL) { fprintf(stderr, "%d: transmitted particle %d is already here...\n", this_node, p->p.identity); errexit(); } } update_local_particles(particles); PART_TRACE(fprintf(stderr, "%d: recv_particles expecting %d bond ints\n", this_node, local_dyn.n)); if (local_dyn.n > 0) { alloc_intlist(&local_dyn, local_dyn.n); MPI_Recv(local_dyn.e, local_dyn.n*sizeof(int), MPI_BYTE, node, REQ_SNDRCV_PART, comm_cart, MPI_STATUS_IGNORE); } read = 0; for (pc = particles->n - transfer; pc < particles->n; pc++) { Particle *p = &particles->part[pc]; if (p->bl.n > 0) { alloc_intlist(&p->bl, p->bl.n); memcpy(p->bl.e, &local_dyn.e[read], p->bl.n*sizeof(int)); read += p->bl.n; } else p->bl.e = NULL; #ifdef EXCLUSIONS if (p->el.n > 0) { alloc_intlist(&p->el, p->el.n); memcpy(p->el.e, &local_dyn.e[read], p->el.n*sizeof(int)); read += p->el.n; } else p->el.e = NULL; #endif } if (local_dyn.n > 0) realloc_intlist(&local_dyn, 0); }