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 203 204 205 206 207 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 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722
|
/* ----------------------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
http://lammps.sandia.gov, Sandia National Laboratories
Steve Plimpton, sjplimp@sandia.gov
Copyright (2003) Sandia Corporation. Under the terms of Contract
DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
certain rights in this software. This software is distributed under
the GNU General Public License.
See the README file in the top-level LAMMPS directory.
------------------------------------------------------------------------- */
#include "lmptype.h"
#include "mpi.h"
#include "stdlib.h"
#include "string.h"
#include "irregular.h"
#include "atom.h"
#include "atom_vec.h"
#include "domain.h"
#include "comm.h"
#include "memory.h"
using namespace LAMMPS_NS;
#define BUFFACTOR 1.5
#define BUFMIN 1000
#define BUFEXTRA 1000
/* ---------------------------------------------------------------------- */
Irregular::Irregular(LAMMPS *lmp) : Pointers(lmp)
{
MPI_Comm_rank(world,&me);
MPI_Comm_size(world,&nprocs);
triclinic = domain->triclinic;
map_style = atom->map_style;
procgrid = comm->procgrid;
grid2proc = comm->grid2proc;
aplan = NULL;
dplan = NULL;
// initialize buffers for atom comm, not used for datum comm
// these can persist for multiple irregular operations
maxsend = BUFMIN;
memory->create(buf_send,maxsend+BUFEXTRA,"comm:buf_send");
maxrecv = BUFMIN;
memory->create(buf_recv,maxrecv,"comm:buf_recv");
}
/* ---------------------------------------------------------------------- */
Irregular::~Irregular()
{
if (aplan) destroy_atom();
if (dplan) destroy_data();
memory->destroy(buf_send);
memory->destroy(buf_recv);
}
/* ----------------------------------------------------------------------
communicate atoms to new owning procs via irregular communication
can be used in place of comm->exchange()
unlike exchange(), allows atoms to have moved arbitrarily long distances
sets up irregular plan, invokes it, destroys it
atoms must be remapped to be inside simulation box before this is called
for triclinic: atoms must be in lamda coords (0-1) before this is called
------------------------------------------------------------------------- */
void Irregular::migrate_atoms()
{
// clear global->local map since atoms move to new procs
// clear old ghosts so map_set() at end will operate only on local atoms
// exchange() doesn't need to clear ghosts b/c borders()
// is called right after and it clears ghosts and calls map_set()
if (map_style) atom->map_clear();
atom->nghost = 0;
atom->avec->clear_bonus();
// subbox bounds for orthogonal or triclinic box
// other comm/domain data used by coord2proc()
double *sublo,*subhi;
if (triclinic == 0) {
sublo = domain->sublo;
subhi = domain->subhi;
} else {
sublo = domain->sublo_lamda;
subhi = domain->subhi_lamda;
}
uniform = comm->uniform;
xsplit = comm->xsplit;
ysplit = comm->ysplit;
zsplit = comm->zsplit;
boxlo = domain->boxlo;
prd = domain->prd;
// loop over atoms, flag any that are not in my sub-box
// fill buffer with atoms leaving my box, using < and >=
// assign which proc it belongs to via coord2proc()
// if coord2proc() returns me, due to round-off
// in triclinic x2lamda(), then keep atom and don't send
// when atom is deleted, fill it in with last atom
AtomVec *avec = atom->avec;
double **x = atom->x;
int nlocal = atom->nlocal;
int nsend = 0;
int nsendatom = 0;
int *sizes = new int[nlocal];
int *proclist = new int[nlocal];
int i = 0;
while (i < nlocal) {
if (x[i][0] < sublo[0] || x[i][0] >= subhi[0] ||
x[i][1] < sublo[1] || x[i][1] >= subhi[1] ||
x[i][2] < sublo[2] || x[i][2] >= subhi[2]) {
proclist[nsendatom] = coord2proc(x[i]);
if (proclist[nsendatom] != me) {
if (nsend > maxsend) grow_send(nsend,1);
sizes[nsendatom] = avec->pack_exchange(i,&buf_send[nsend]);
nsend += sizes[nsendatom];
nsendatom++;
avec->copy(nlocal-1,i,1);
nlocal--;
} else i++;
} else i++;
}
atom->nlocal = nlocal;
// create irregular communication plan, perform comm, destroy plan
// returned nrecv = size of buffer needed for incoming atoms
int nrecv = create_atom(nsendatom,sizes,proclist);
if (nrecv > maxrecv) grow_recv(nrecv);
exchange_atom(buf_send,sizes,buf_recv);
destroy_atom();
delete [] sizes;
delete [] proclist;
// add received atoms to my list
int m = 0;
while (m < nrecv) m += avec->unpack_exchange(&buf_recv[m]);
// reset global->local map
if (map_style) atom->map_set();
}
/* ----------------------------------------------------------------------
create a communication plan for atoms
n = # of atoms to send
sizes = # of doubles for each atom
proclist = proc to send each atom to (not including self)
return total # of doubles I will recv (not including self)
------------------------------------------------------------------------- */
int Irregular::create_atom(int n, int *sizes, int *proclist)
{
int i;
// allocate plan and work vectors
if (aplan) destroy_atom();
aplan = (PlanAtom *) memory->smalloc(sizeof(PlanAtom),"irregular:aplan");
int *list = new int[nprocs];
int *count = new int[nprocs];
// nrecv = # of messages I receive
for (i = 0; i < nprocs; i++) {
list[i] = 0;
count[i] = 1;
}
for (i = 0; i < n; i++) list[proclist[i]] = 1;
int nrecv;
MPI_Reduce_scatter(list,&nrecv,count,MPI_INT,MPI_SUM,world);
// allocate receive arrays
int *proc_recv = new int[nrecv];
int *length_recv = new int[nrecv];
MPI_Request *request = new MPI_Request[nrecv];
MPI_Status *status = new MPI_Status[nrecv];
// nsend = # of messages I send
for (i = 0; i < nprocs; i++) list[i] = 0;
for (i = 0; i < n; i++) list[proclist[i]] += sizes[i];
int nsend = 0;
for (i = 0; i < nprocs; i++)
if (list[i]) nsend++;
// allocate send arrays
int *proc_send = new int[nsend];
int *length_send = new int[nsend];
int *num_send = new int[nsend];
int *index_send = new int[n];
int *offset_send = new int[n];
// list still stores size of message for procs I send to
// proc_send = procs I send to
// length_send = total size of message I send to each proc
// to balance pattern of send messages:
// each proc begins with iproc > me, continues until iproc = me
// reset list to store which send message each proc corresponds to
int iproc = me;
int isend = 0;
for (i = 0; i < nprocs; i++) {
iproc++;
if (iproc == nprocs) iproc = 0;
if (list[iproc] > 0) {
proc_send[isend] = iproc;
length_send[isend] = list[iproc];
list[iproc] = isend;
isend++;
}
}
// num_send = # of atoms I send to each proc
for (i = 0; i < nsend; i++) num_send[i] = 0;
for (i = 0; i < n; i++) {
isend = list[proclist[i]];
num_send[isend]++;
}
// count = offsets into index_send for each proc I send to
// index_send = list of which atoms to send to each proc
// 1st N1 values are atom indices for 1st proc,
// next N2 values are atom indices for 2nd proc, etc
// offset_send = where each atom starts in send buffer
count[0] = 0;
for (i = 1; i < nsend; i++) count[i] = count[i-1] + num_send[i-1];
for (i = 0; i < n; i++) {
isend = list[proclist[i]];
index_send[count[isend]++] = i;
if (i) offset_send[i] = offset_send[i-1] + sizes[i-1];
else offset_send[i] = 0;
}
// tell receivers how much data I send
// sendmax = largest # of doubles I send in a single message
int sendmax = 0;
for (i = 0; i < nsend; i++) {
MPI_Send(&length_send[i],1,MPI_INT,proc_send[i],0,world);
sendmax = MAX(sendmax,length_send[i]);
}
// receive incoming messages
// proc_recv = procs I recv from
// length_recv = total size of message each proc sends me
// nrecvsize = total size of data I recv
int nrecvsize = 0;
for (i = 0; i < nrecv; i++) {
MPI_Recv(&length_recv[i],1,MPI_INT,MPI_ANY_SOURCE,0,world,status);
proc_recv[i] = status->MPI_SOURCE;
nrecvsize += length_recv[i];
}
// barrier to insure all MPI_ANY_SOURCE messages are received
// else another proc could proceed to exchange_atom() and send to me
MPI_Barrier(world);
// free work vectors
delete [] count;
delete [] list;
// initialize plan
aplan->nsend = nsend;
aplan->nrecv = nrecv;
aplan->sendmax = sendmax;
aplan->proc_send = proc_send;
aplan->length_send = length_send;
aplan->num_send = num_send;
aplan->index_send = index_send;
aplan->offset_send = offset_send;
aplan->proc_recv = proc_recv;
aplan->length_recv = length_recv;
aplan->request = request;
aplan->status = status;
return nrecvsize;
}
/* ----------------------------------------------------------------------
communicate atoms via PlanAtom
sendbuf = list of atoms to send
sizes = # of doubles for each atom
recvbuf = received atoms
------------------------------------------------------------------------- */
void Irregular::exchange_atom(double *sendbuf, int *sizes, double *recvbuf)
{
int i,m,n,offset,num_send;
// post all receives
offset = 0;
for (int irecv = 0; irecv < aplan->nrecv; irecv++) {
MPI_Irecv(&recvbuf[offset],aplan->length_recv[irecv],MPI_DOUBLE,
aplan->proc_recv[irecv],0,world,&aplan->request[irecv]);
offset += aplan->length_recv[irecv];
}
// allocate buf for largest send
double *buf;
memory->create(buf,aplan->sendmax,"irregular:buf");
// send each message
// pack buf with list of atoms
// m = index of atom in sendbuf
int *index_send = aplan->index_send;
int nsend = aplan->nsend;
n = 0;
for (int isend = 0; isend < nsend; isend++) {
offset = 0;
num_send = aplan->num_send[isend];
for (i = 0; i < num_send; i++) {
m = index_send[n++];
memcpy(&buf[offset],&sendbuf[aplan->offset_send[m]],
sizes[m]*sizeof(double));
offset += sizes[m];
}
MPI_Send(buf,aplan->length_send[isend],MPI_DOUBLE,
aplan->proc_send[isend],0,world);
}
// free temporary send buffer
memory->destroy(buf);
// wait on all incoming messages
if (aplan->nrecv) MPI_Waitall(aplan->nrecv,aplan->request,aplan->status);
}
/* ----------------------------------------------------------------------
destroy communication plan for atoms
------------------------------------------------------------------------- */
void Irregular::destroy_atom()
{
delete [] aplan->proc_send;
delete [] aplan->length_send;
delete [] aplan->num_send;
delete [] aplan->index_send;
delete [] aplan->offset_send;
delete [] aplan->proc_recv;
delete [] aplan->length_recv;
delete [] aplan->request;
delete [] aplan->status;
memory->sfree(aplan);
aplan = NULL;
}
/* ----------------------------------------------------------------------
create a communication plan for datums
n = # of datums to send
proclist = proc to send each datum to (including self)
return total # of datums I will recv (including self)
------------------------------------------------------------------------- */
int Irregular::create_data(int n, int *proclist)
{
int i,m;
// allocate plan and work vectors
dplan = (PlanData *) memory->smalloc(sizeof(PlanData),"irregular:dplan");
int *list = new int[nprocs];
int *count = new int[nprocs];
// nrecv = # of messages I receive
for (i = 0; i < nprocs; i++) {
list[i] = 0;
count[i] = 1;
}
for (i = 0; i < n; i++) list[proclist[i]] = 1;
int nrecv;
MPI_Reduce_scatter(list,&nrecv,count,MPI_INT,MPI_SUM,world);
if (list[me]) nrecv--;
// allocate receive arrays
int *proc_recv = new int[nrecv];
int *num_recv = new int[nrecv];
MPI_Request *request = new MPI_Request[nrecv];
MPI_Status *status = new MPI_Status[nrecv];
// nsend = # of messages I send
for (i = 0; i < nprocs; i++) list[i] = 0;
for (i = 0; i < n; i++) list[proclist[i]]++;
int nsend = 0;
for (i = 0; i < nprocs; i++)
if (list[i]) nsend++;
if (list[me]) nsend--;
// allocate send and self arrays
int *proc_send = new int[nsend];
int *num_send = new int[nsend];
int *index_send = new int[n-list[me]];
int *index_self = new int[list[me]];
// proc_send = procs I send to
// num_send = # of datums I send to each proc
// num_self = # of datums I copy to self
// to balance pattern of send messages:
// each proc begins with iproc > me, continues until iproc = me
// reset list to store which send message each proc corresponds to
int num_self;
int iproc = me;
int isend = 0;
for (i = 0; i < nprocs; i++) {
iproc++;
if (iproc == nprocs) iproc = 0;
if (iproc == me) num_self = list[iproc];
else if (list[iproc] > 0) {
proc_send[isend] = iproc;
num_send[isend] = list[iproc];
list[iproc] = isend;
isend++;
}
}
list[me] = 0;
// count = offsets into index_send for each proc I send to
// m = ptr into index_self
// index_send = list of which datums to send to each proc
// 1st N1 values are datum indices for 1st proc,
// next N2 values are datum indices for 2nd proc, etc
count[0] = 0;
for (i = 1; i < nsend; i++) count[i] = count[i-1] + num_send[i-1];
m = 0;
for (i = 0; i < n; i++) {
iproc = proclist[i];
if (iproc == me) index_self[m++] = i;
else {
isend = list[iproc];
index_send[count[isend]++] = i;
}
}
// tell receivers how much data I send
// sendmax = largest # of datums I send in a single message
int sendmax = 0;
for (i = 0; i < nsend; i++) {
MPI_Send(&num_send[i],1,MPI_INT,proc_send[i],0,world);
sendmax = MAX(sendmax,num_send[i]);
}
// receive incoming messages
// proc_recv = procs I recv from
// num_recv = total size of message each proc sends me
// nrecvsize = total size of data I recv
int nrecvsize = 0;
for (i = 0; i < nrecv; i++) {
MPI_Recv(&num_recv[i],1,MPI_INT,MPI_ANY_SOURCE,0,world,status);
proc_recv[i] = status->MPI_SOURCE;
nrecvsize += num_recv[i];
}
nrecvsize += num_self;
// barrier to insure all MPI_ANY_SOURCE messages are received
// else another proc could proceed to exchange_data() and send to me
MPI_Barrier(world);
// free work vectors
delete [] count;
delete [] list;
// initialize plan and return it
dplan->nsend = nsend;
dplan->nrecv = nrecv;
dplan->sendmax = sendmax;
dplan->proc_send = proc_send;
dplan->num_send = num_send;
dplan->index_send = index_send;
dplan->proc_recv = proc_recv;
dplan->num_recv = num_recv;
dplan->num_self = num_self;
dplan->index_self = index_self;
dplan->request = request;
dplan->status = status;
return nrecvsize;
}
/* ----------------------------------------------------------------------
communicate datums via PlanData
sendbuf = list of datums to send
nbytes = size of each datum
recvbuf = received datums (including copied from me)
------------------------------------------------------------------------- */
void Irregular::exchange_data(char *sendbuf, int nbytes, char *recvbuf)
{
int i,m,n,offset,num_send;
// post all receives, starting after self copies
offset = dplan->num_self*nbytes;
for (int irecv = 0; irecv < dplan->nrecv; irecv++) {
MPI_Irecv(&recvbuf[offset],dplan->num_recv[irecv]*nbytes,MPI_CHAR,
dplan->proc_recv[irecv],0,world,&dplan->request[irecv]);
offset += dplan->num_recv[irecv]*nbytes;
}
// allocate buf for largest send
char *buf;
memory->create(buf,dplan->sendmax*nbytes,"irregular:buf");
// send each message
// pack buf with list of datums
// m = index of datum in sendbuf
int *index_send = dplan->index_send;
int nsend = dplan->nsend;
n = 0;
for (int isend = 0; isend < nsend; isend++) {
num_send = dplan->num_send[isend];
for (i = 0; i < num_send; i++) {
m = index_send[n++];
memcpy(&buf[i*nbytes],&sendbuf[m*nbytes],nbytes);
}
MPI_Send(buf,dplan->num_send[isend]*nbytes,MPI_CHAR,
dplan->proc_send[isend],0,world);
}
// free temporary send buffer
memory->destroy(buf);
// copy datums to self, put at beginning of recvbuf
int *index_self = dplan->index_self;
int num_self = dplan->num_self;
for (i = 0; i < num_self; i++) {
m = index_self[i];
memcpy(&recvbuf[i*nbytes],&sendbuf[m*nbytes],nbytes);
}
// wait on all incoming messages
if (dplan->nrecv) MPI_Waitall(dplan->nrecv,dplan->request,dplan->status);
}
/* ----------------------------------------------------------------------
destroy communication plan for datums
------------------------------------------------------------------------- */
void Irregular::destroy_data()
{
delete [] dplan->proc_send;
delete [] dplan->num_send;
delete [] dplan->index_send;
delete [] dplan->proc_recv;
delete [] dplan->num_recv;
delete [] dplan->index_self;
delete [] dplan->request;
delete [] dplan->status;
memory->sfree(dplan);
dplan = NULL;
}
/* ----------------------------------------------------------------------
determine which proc owns atom with coord x[3]
x will be in box (orthogonal) or lamda coords (triclinic)
for uniform = 1, directly calculate owning proc
for non-uniform, iteratively find owning proc via binary search
------------------------------------------------------------------------- */
int Irregular::coord2proc(double *x)
{
int loc[3];
if (uniform) {
if (triclinic == 0) {
loc[0] = static_cast<int> (procgrid[0] * (x[0]-boxlo[0]) / prd[0]);
loc[1] = static_cast<int> (procgrid[1] * (x[1]-boxlo[1]) / prd[1]);
loc[2] = static_cast<int> (procgrid[2] * (x[2]-boxlo[2]) / prd[2]);
} else {
loc[0] = static_cast<int> (procgrid[0] * x[0]);
loc[1] = static_cast<int> (procgrid[1] * x[1]);
loc[2] = static_cast<int> (procgrid[2] * x[2]);
}
} else {
if (triclinic == 0) {
loc[0] = binary((x[0]-boxlo[0])/prd[0],procgrid[0],xsplit);
loc[1] = binary((x[1]-boxlo[1])/prd[1],procgrid[1],ysplit);
loc[2] = binary((x[2]-boxlo[2])/prd[2],procgrid[2],zsplit);
} else {
loc[0] = binary(x[0],procgrid[0],xsplit);
loc[1] = binary(x[1],procgrid[1],ysplit);
loc[2] = binary(x[2],procgrid[2],zsplit);
}
}
if (loc[0] < 0) loc[0] = 0;
if (loc[0] >= procgrid[0]) loc[0] = procgrid[0] - 1;
if (loc[1] < 0) loc[1] = 0;
if (loc[1] >= procgrid[1]) loc[1] = procgrid[1] - 1;
if (loc[2] < 0) loc[2] = 0;
if (loc[2] >= procgrid[2]) loc[2] = procgrid[2] - 1;
return grid2proc[loc[0]][loc[1]][loc[2]];
}
/* ----------------------------------------------------------------------
binary search for value in N-length ascending vec
value may be outside range of vec limits
always return index from 0 to N-1 inclusive
return 0 if value < vec[0]
reutrn N-1 if value >= vec[N-1]
return index = 1 to N-2 if vec[index] <= value < vec[index+1]
------------------------------------------------------------------------- */
int Irregular::binary(double value, int n, double *vec)
{
int lo = 0;
int hi = n-1;
if (value < vec[lo]) return lo;
if (value >= vec[hi]) return hi;
// insure vec[lo] <= value < vec[hi] at every iteration
// done when lo,hi are adjacent
int index = (lo+hi)/2;
while (lo < hi-1) {
if (value < vec[index]) hi = index;
else if (value >= vec[index]) lo = index;
index = (lo+hi)/2;
}
return index;
}
/* ----------------------------------------------------------------------
realloc the size of the send buffer as needed with BUFFACTOR & BUFEXTRA
if flag = 1, realloc
if flag = 0, don't need to realloc with copy, just free/malloc
------------------------------------------------------------------------- */
void Irregular::grow_send(int n, int flag)
{
maxsend = static_cast<int> (BUFFACTOR * n);
if (flag)
memory->grow(buf_send,maxsend+BUFEXTRA,"comm:buf_send");
else {
memory->destroy(buf_send);
memory->create(buf_send,maxsend+BUFEXTRA,"comm:buf_send");
}
}
/* ----------------------------------------------------------------------
free/malloc the size of the recv buffer as needed with BUFFACTOR
------------------------------------------------------------------------- */
void Irregular::grow_recv(int n)
{
maxrecv = static_cast<int> (BUFFACTOR * n);
memory->destroy(buf_recv);
memory->create(buf_recv,maxrecv,"comm:buf_recv");
}
/* ----------------------------------------------------------------------
return # of bytes of allocated memory
------------------------------------------------------------------------- */
bigint Irregular::memory_usage()
{
bigint bytes = memory->usage(buf_send,maxsend);
bytes += memory->usage(buf_recv,maxrecv);
return bytes;
}
|