git-svn-id: svn://svn.icms.temple.edu/lammps-ro/trunk@12239 f3b2605a-c512-4ea7-a41b-209d697bcdaa

This commit is contained in:
sjplimp 2014-08-01 19:15:50 +00:00
parent df5c028ad4
commit 5ac9e9a785
6 changed files with 423 additions and 58 deletions

View File

@ -488,10 +488,14 @@ int *Balance::bisection(int sortflag)
rcb->invert(sortflag);
// NOTE: this logic is specific to orthogonal boxes, not triclinic
// insure comm->rcbcut will be exactly equal to sub-box boundaries it touches
comm->rcbnew = 1;
comm->rcbcut = rcb->cut;
comm->rcbcutdim = rcb->cutdim;
int idim = rcb->cutdim;
double split = (rcb->cut - boxlo[idim]) / prd[idim];
comm->rcbcut = boxlo[idim] + prd[idim]*split;
comm->rcbcutdim = idim;
double (*mysplit)[2] = comm->mysplit;

View File

@ -678,12 +678,14 @@ void CommBrick::exchange()
atom->nlocal = nlocal;
// send/recv atoms in both directions
// send size of message first so receiver can realloc buf_recv if needed
// if 1 proc in dimension, no send/recv, set recv buf to send buf
// if 2 procs in dimension, single send/recv
// if more than 2 procs in dimension, send/recv to both neighbors
if (procgrid[dim] == 1) {
nrecv = nsend;
nrecv = nsend; // NOTE: could just be nrecv = 0 ??
// no buf, just buf_recv
buf = buf_send;
} else {
@ -714,6 +716,7 @@ void CommBrick::exchange()
// check incoming atoms to see if they are in my box
// if so, add to my list
// check is only for this dimension, may be passed to another proc
m = 0;
while (m < nrecv) {

View File

@ -74,6 +74,7 @@ CommTiled::~CommTiled()
memory->destroy(buf_recv);
memory->destroy(overlap);
deallocate_swap(nswap);
memory->sfree(rcbinfo);
}
/* ----------------------------------------------------------------------
@ -96,6 +97,8 @@ void CommTiled::init_buffers()
nswap = 2 * domain->dimension;
allocate_swap(nswap);
rcbinfo = NULL;
}
/* ----------------------------------------------------------------------
@ -173,9 +176,9 @@ void CommTiled::init()
void CommTiled::setup()
{
int i;
int i,n;
// domain properties used in setup and methods it calls
// domain properties used in setup method and methods it calls
prd = domain->prd;
boxlo = domain->boxlo;
@ -191,19 +194,25 @@ void CommTiled::setup()
if (layout != LAYOUT_TILED) {
box_drop = &CommTiled::box_drop_brick;
box_other = &CommTiled::box_other_brick;
box_touch = &CommTiled::box_touch_brick;
} else {
box_drop = &CommTiled::box_drop_tiled;
box_other = &CommTiled::box_other_tiled;
box_touch = &CommTiled::box_touch_tiled;
}
// if RCB decomp has just changed, gather needed global RCB info
// if RCB decomp exists and just changed, gather needed global RCB info
if (rcbnew) {
if (!rcbinfo)
rcbinfo = (RCBinfo *)
memory->smalloc(nprocs*sizeof(RCBinfo),"comm:rcbinfo");
rcbnew = 0;
memcpy(&rcbinfo[me].mysplit[0][0],&mysplit[0][0],6*sizeof(double));
rcbinfo[me].cut = rcbcut;
rcbinfo[me].dim = rcbcutdim;
MPI_Allgather(&rcbinfo[me],sizeof(RCBinfo),MPI_CHAR,
RCBinfo rcbone;
memcpy(&rcbone.mysplit[0][0],&mysplit[0][0],6*sizeof(double));
rcbone.cut = rcbcut;
rcbone.dim = rcbcutdim;
MPI_Allgather(&rcbone,sizeof(RCBinfo),MPI_CHAR,
rcbinfo,sizeof(RCBinfo),MPI_CHAR,world);
}
@ -218,11 +227,13 @@ void CommTiled::setup()
error->all(FLERR,"Communication cutoff for comm_style tiled "
"cannot exceed periodic box length");
// setup forward/reverse communication
// loop over 6 swap directions
// determine which procs I will send to and receive from in each swap
// done by intersecting ghost box with all proc sub-boxes it overlaps
// sets nsendproc, nrecvproc, sendproc, recvproc
// sets sendother, sendself, pbc_flag, pbc, sendbox
// resets nprocmax
int noverlap1,indexme;
double lo1[3],hi1[3],lo2[3],hi2[3];
@ -292,12 +303,12 @@ void CommTiled::setup()
// reallocate 2nd dimensions of all send/recv arrays, based on noverlap
// # of sends of this swap = # of recvs of nswap +/- 1
if (noverlap > nprocmax[iswap]) {
int oldmax = nprocmax[iswap];
while (nprocmax[iswap] < noverlap) nprocmax[iswap] += DELTA_PROCS;
grow_swap_send(iswap,nprocmax[iswap],oldmax);
if (iswap == 0) grow_swap_recv(iswap+1,nprocmax[iswap]);
else grow_swap_recv(iswap-1,nprocmax[iswap]);
if (noverlap > nprocmax[nswap]) {
int oldmax = nprocmax[nswap];
while (nprocmax[nswap] < noverlap) nprocmax[nswap] += DELTA_PROCS;
grow_swap_send(nswap,nprocmax[nswap],oldmax);
if (iswap == 0) grow_swap_recv(nswap+1,nprocmax[nswap]);
else grow_swap_recv(nswap-1,nprocmax[nswap]);
}
// overlap how has list of noverlap procs
@ -308,6 +319,11 @@ void CommTiled::setup()
if (noverlap-sendself[nswap]) sendother[nswap] = 1;
else sendother[nswap] = 0;
//MPI_Barrier(world);
//printf("AAA nswap %d me %d: noverlap %d: %g %g: %g %g\n",
// nswap,me,noverlap,sublo[0],sublo[1],subhi[0],subhi[1]);
//if (nswap == 0) error->all(FLERR,"ALL DONE");
nsendproc[nswap] = noverlap;
for (i = 0; i < noverlap; i++) sendproc[nswap][i] = overlap[i];
if (iswap == 0) {
@ -357,19 +373,21 @@ void CommTiled::setup()
if (iswap == 0) {
sbox[idim] = sublo[idim];
sbox[3+idim] = MIN(sublo[idim]+cut,subhi[idim]);
if (i < noverlap1) sbox[3+idim] = MIN(sbox[3+idim]+cut,subhi[idim]);
else sbox[3+idim] = MIN(sbox[3+idim]-prd[idim]+cut,subhi[idim]);
} else {
sbox[idim] = MAX(subhi[idim]-cut,sublo[idim]);
if (i < noverlap1) sbox[idim] = MAX(sbox[idim]-cut,sublo[idim]);
else sbox[idim] = MAX(sbox[idim]+prd[idim]-cut,sublo[idim]);
sbox[3+idim] = subhi[idim];
}
if (idim >= 1) {
if (sbox[0] == sublo[0]) sbox[0] -= cut;
if (sbox[4] == subhi[0]) sbox[4] += cut;
if (sbox[3] == subhi[0]) sbox[3] += cut;
}
if (idim == 2) {
if (sbox[1] == sublo[1]) sbox[1] -= cut;
if (sbox[5] == subhi[1]) sbox[5] += cut;
if (sbox[4] == subhi[1]) sbox[4] += cut;
}
memcpy(sendbox[nswap][i],sbox,6*sizeof(double));
@ -378,6 +396,112 @@ void CommTiled::setup()
nswap++;
}
}
// setup exchange communication = subset of forward/reverse comm
// loop over 6 swap directions
// determine which procs I will send to and receive from in each swap
// subset of procs that touch my proc in forward/reverse comm
// sets nesendproc, nerecvproc, esendproc, erecvproc
// resets neprocmax
// NOTE: should there be a unique neprocmax?
printf("SUBBOX %d: %g %g: %g %g\n",me,sublo[0],sublo[1],subhi[0],subhi[1]);
MPI_Barrier(world);
nswap = 0;
for (int idim = 0; idim < dimension; idim++) {
for (int iswap = 0; iswap < 2; iswap++) {
noverlap = 0;
n = nsendproc[nswap];
for (i = 0; i < n; i++) {
if (sendproc[nswap][i] == me) continue;
if ((this->*box_touch)(sendproc[nswap][i],idim,iswap))
overlap[noverlap++] = sendproc[nswap][i];
}
// reallocate esendproc or erecvproc if needed based on novlerap
if (noverlap > neprocmax[nswap]) {
while (neprocmax[nswap] < noverlap) neprocmax[nswap] += DELTA_PROCS;
n = neprocmax[nswap];
delete [] esendproc[nswap];
esendproc[nswap] = new int[n];
if (iswap == 0) {
delete [] erecvproc[nswap+1];
erecvproc[nswap+1] = new int[n];
} else {
delete [] erecvproc[nswap-1];
erecvproc[nswap-1] = new int[n];
}
}
nesendproc[nswap] = noverlap;
for (i = 0; i < noverlap; i++) esendproc[nswap][i] = overlap[i];
if (iswap == 0) {
nerecvproc[nswap+1] = noverlap;
for (i = 0; i < noverlap; i++) erecvproc[nswap+1][i] = overlap[i];
} else {
nerecvproc[nswap-1] = noverlap;
for (i = 0; i < noverlap; i++) erecvproc[nswap-1][i] = overlap[i];
}
nswap++;
}
}
// reallocate MPI Requests and Statuses as needed
int nmax = 0;
for (i = 0; i < nswap; i++) nmax = MAX(nmax,nprocmax[i]);
if (nmax > maxreqstat) {
maxreqstat = nmax;
delete [] requests;
delete [] statuses;
requests = new MPI_Request[maxreqstat];
statuses = new MPI_Status[maxreqstat];
}
// DEBUG
/*
for (i = 0; i < nswap; i++) {
if (nsendproc[i] == 1)
printf("SETUP SEND %d %d: nsend %d self %d sproc0 %d: "
"%g %g %g: %g %g %g\n",
i,me,nsendproc[i],sendself[i],sendproc[i][0],
sendbox[i][0][0],
sendbox[i][0][1],
sendbox[i][0][2],
sendbox[i][0][3],
sendbox[i][0][4],
sendbox[i][0][5]);
else
printf("SETUP SEND %d %d: nsend %d self %d sprocs %d %d: "
"%g %g %g: %g %g %g: %g %g %g: %g %g %g\n",
i,me,nsendproc[i],sendself[i],sendproc[i][0],sendproc[i][1],
sendbox[i][0][0],
sendbox[i][0][1],
sendbox[i][0][2],
sendbox[i][0][3],
sendbox[i][0][4],
sendbox[i][0][5],
sendbox[i][1][0],
sendbox[i][1][1],
sendbox[i][1][2],
sendbox[i][1][3],
sendbox[i][1][4],
sendbox[i][1][5]);
if (nrecvproc[i] == 1)
printf("SETUP RECV %d %d: nrecv %d other %d rproc0 %d\n",
i,me,nrecvproc[i],sendother[i],recvproc[i][0]);
else
printf("SETUP RECV %d %d: nrecv %d other %d rprocs %d %d\n",
i,me,nrecvproc[i],sendother[i],recvproc[i][0],recvproc[i][1]);
}
MPI_Barrier(world);
*/
}
/* ----------------------------------------------------------------------
@ -426,7 +550,7 @@ void CommTiled::forward_comm(int dummy)
} else if (ghost_velocity) {
if (sendother[iswap]) {
for (i = 0; i < nrecv; i++)
MPI_Irecv(&buf_recv[forward_recv_offset[iswap][i]],
MPI_Irecv(&buf_recv[size_forward*forward_recv_offset[iswap][i]],
size_forward_recv[iswap][i],
MPI_DOUBLE,recvproc[iswap][i],0,world,&requests[i]);
for (i = 0; i < nsend; i++) {
@ -447,14 +571,15 @@ void CommTiled::forward_comm(int dummy)
for (i = 0; i < nrecv; i++) {
MPI_Waitany(nrecv,requests,&irecv,&status);
avec->unpack_comm_vel(recvnum[iswap][irecv],firstrecv[iswap][irecv],
&buf_recv[forward_recv_offset[iswap][irecv]]);
&buf_recv[size_forward*
forward_recv_offset[iswap][irecv]]);
}
}
} else {
if (sendother[iswap]) {
for (i = 0; i < nrecv; i++)
MPI_Irecv(&buf_recv[forward_recv_offset[iswap][i]],
MPI_Irecv(&buf_recv[size_forward*forward_recv_offset[iswap][i]],
size_forward_recv[iswap][i],
MPI_DOUBLE,recvproc[iswap][i],0,world,&requests[i]);
for (i = 0; i < nsendproc[iswap]; i++) {
@ -475,7 +600,8 @@ void CommTiled::forward_comm(int dummy)
for (i = 0; i < nrecv; i++) {
MPI_Waitany(nrecv,requests,&irecv,&status);
avec->unpack_comm(recvnum[iswap][irecv],firstrecv[iswap][irecv],
&buf_recv[forward_recv_offset[iswap][irecv]]);
&buf_recv[size_forward*
forward_recv_offset[iswap][irecv]]);
}
}
}
@ -509,7 +635,7 @@ void CommTiled::reverse_comm()
if (comm_f_only) {
if (sendother[iswap]) {
for (i = 0; i < nsend; i++)
MPI_Irecv(&buf_recv[reverse_recv_offset[iswap][i]],
MPI_Irecv(&buf_recv[size_reverse*reverse_recv_offset[iswap][i]],
size_reverse_recv[iswap][i],MPI_DOUBLE,
sendproc[iswap][i],0,world,&requests[i]);
for (i = 0; i < nrecv; i++)
@ -526,14 +652,15 @@ void CommTiled::reverse_comm()
for (i = 0; i < nsend; i++) {
MPI_Waitany(nsend,requests,&irecv,&status);
avec->unpack_reverse(sendnum[iswap][irecv],sendlist[iswap][irecv],
&buf_recv[reverse_recv_offset[iswap][irecv]]);
&buf_recv[size_reverse*
reverse_recv_offset[iswap][irecv]]);
}
}
} else {
if (sendother[iswap]) {
for (i = 0; i < nsend; i++)
MPI_Irecv(&buf_recv[reverse_recv_offset[iswap][i]],
MPI_Irecv(&buf_recv[size_reverse*reverse_recv_offset[iswap][i]],
size_reverse_recv[iswap][i],MPI_DOUBLE,
sendproc[iswap][i],0,world,&requests[i]);
for (i = 0; i < nrecv; i++) {
@ -554,7 +681,8 @@ void CommTiled::reverse_comm()
for (i = 0; i < nsend; i++) {
MPI_Waitany(nsend,requests,&irecv,&status);
avec->unpack_reverse(sendnum[iswap][irecv],sendlist[iswap][irecv],
&buf_recv[reverse_recv_offset[iswap][irecv]]);
&buf_recv[size_reverse*
reverse_recv_offset[iswap][irecv]]);
}
}
}
@ -563,6 +691,7 @@ void CommTiled::reverse_comm()
/* ----------------------------------------------------------------------
exchange: move atoms to correct processors
NOTE: need to re-doc this
atoms exchanged with all 6 stencil neighbors
send out atoms that have left my box, receive ones entering my box
atoms will be lost if not inside some proc's box
@ -574,14 +703,134 @@ void CommTiled::reverse_comm()
void CommTiled::exchange()
{
// loop over atoms
// if not outside my box, continue
// find which proc it is in
// find which one of my touching procs it is, else lost
// make sure all atoms are "lost" that should be (e.g. outside non-PBC)
// add to list to send to that proc
// loop over touching procs
// send buffer to them
int i,m,nsend,nrecv,nsendsize,nrecvsize,nlocal,dim,offset;
double lo,hi,value;
double **x;
AtomVec *avec = atom->avec;
MPI_Barrier(world);
printf("PREEXCH %d %d\n",me,atom->nlocal);
MPI_Barrier(world);
// clear global->local map for owned and ghost atoms
// b/c atoms migrate to new procs in exchange() and
// new ghosts are created in borders()
// map_set() is done at end of borders()
// clear ghost count and any ghost bonus data internal to AtomVec
if (map_style) atom->map_clear();
atom->nghost = 0;
atom->avec->clear_bonus();
// insure send buf is large enough for single atom
// fixes can change per-atom size requirement on-the-fly
int bufextra_old = bufextra;
maxexchange = maxexchange_atom + maxexchange_fix;
bufextra = maxexchange + BUFEXTRA;
if (bufextra > bufextra_old)
memory->grow(buf_send,maxsend+bufextra,"comm:buf_send");
// subbox bounds for orthogonal or triclinic
if (triclinic == 0) {
sublo = domain->sublo;
subhi = domain->subhi;
} else {
sublo = domain->sublo_lamda;
subhi = domain->subhi_lamda;
}
// loop over swaps in all dimensions
for (int iswap = 0; iswap < nswap; iswap++) {
dim = iswap/2;
// fill buffer with atoms leaving my box, using < and >=
// when atom is deleted, fill it in with last atom
x = atom->x;
lo = sublo[dim];
hi = subhi[dim];
nlocal = atom->nlocal;
i = nsendsize = 0;
if (iswap % 2 == 0) {
while (i < nlocal) {
if (x[i][dim] < lo) {
printf("SEND1 from me %d on swap %d: %d: %24.18g %24.18g\n",
me,iswap,atom->tag[i],x[i][dim],lo);
if (nsendsize > maxsend) grow_send(nsendsize,1);
nsendsize += avec->pack_exchange(i,&buf_send[nsendsize]);
avec->copy(nlocal-1,i,1);
nlocal--;
} else i++;
}
} else {
while (i < nlocal) {
if (x[i][dim] >= hi) {
printf("SEND2 from me %d on swap %d: %d: %24.18g %24.18g\n",
me,iswap,atom->tag[i],x[i][dim],hi);
if (nsendsize > maxsend) grow_send(nsendsize,1);
nsendsize += avec->pack_exchange(i,&buf_send[nsendsize]);
avec->copy(nlocal-1,i,1);
nlocal--;
} else i++;
}
}
atom->nlocal = nlocal;
// send and recv atoms from neighbor procs that touch my sub-box in dim
// no send/recv with self
// send size of message first
// receiver may receive multiple messages, can realloc buf_recv if needed
nsend = nesendproc[iswap];
nrecv = nerecvproc[iswap];
for (m = 0; m < nrecv; m++)
MPI_Irecv(&recvnum[iswap][m],1,MPI_INT,
erecvproc[iswap][m],0,world,&requests[m]);
for (m = 0; m < nsend; m++)
MPI_Send(&nsendsize,1,MPI_INT,esendproc[iswap][m],0,world);
MPI_Waitall(nrecv,requests,statuses);
nrecvsize = 0;
for (m = 0; m < nrecv; m++) nrecvsize += recvnum[iswap][m];
if (nrecvsize > maxrecv) grow_recv(nrecvsize);
offset = 0;
for (m = 0; m < nrecv; m++) {
MPI_Irecv(&buf_recv[offset],recvnum[iswap][m],
MPI_DOUBLE,erecvproc[iswap][m],0,world,&requests[m]);
offset += recvnum[iswap][m];
}
for (m = 0; m < nsend; m++)
MPI_Send(buf_send,nsendsize,MPI_DOUBLE,esendproc[iswap][m],0,world);
MPI_Waitall(nrecv,requests,statuses);
// check incoming atoms to see if they are in my box
// if so, add to my list
// check is only for this dimension, may be passed to another proc
m = 0;
while (m < nrecvsize) {
value = buf_recv[m+dim+1];
if (value >= lo && value < hi) {
m += avec->unpack_exchange(&buf_recv[m]);
printf("RECV from me %d on swap %d: %d\n",me,iswap,
atom->tag[atom->nlocal-1]);
}
else m += static_cast<int> (buf_recv[m]);
}
}
MPI_Barrier(world);
printf("POSTEXCH %d %d\n",me,atom->nlocal);
MPI_Barrier(world);
if (atom->firstgroupname) atom->first_reorder();
}
/* ----------------------------------------------------------------------
@ -609,7 +858,7 @@ void CommTiled::borders()
int smax = 0;
int rmax = 0;
// loop over all swaps in all dimensions
// loop over swaps in all dimensions
for (int iswap = 0; iswap < nswap; iswap++) {
@ -632,18 +881,20 @@ void CommTiled::borders()
else nlast = atom->nlocal + atom->nghost;
ncount = 0;
for (i = 0; i < nlocal; i++)
for (i = 0; i < nlocal; i++) {
if (x[i][0] >= xlo && x[i][0] <= xhi &&
x[i][1] >= ylo && x[i][1] <= yhi &&
x[i][2] >= zlo && x[i][2] <= zhi) {
if (ncount == maxsendlist[iswap][m]) grow_list(iswap,i,ncount);
if (ncount == maxsendlist[iswap][m]) grow_list(iswap,m,ncount);
sendlist[iswap][m][ncount++] = i;
}
}
for (i = atom->nlocal; i < nlast; i++)
if (x[i][0] >= xlo && x[i][0] <= xhi &&
x[i][1] >= ylo && x[i][1] <= yhi &&
x[i][2] >= zlo && x[i][2] <= zhi) {
if (ncount == maxsendlist[iswap][m]) grow_list(iswap,i,ncount);
if (ncount == maxsendlist[iswap][m]) grow_list(iswap,m,ncount);
sendlist[iswap][m][ncount++] = i;
}
sendnum[iswap][m] = ncount;
@ -668,12 +919,18 @@ void CommTiled::borders()
// setup other per swap/proc values from sendnum and recvnum
for (m = 0; m < nsendproc[iswap]; m++) {
size_reverse_recv[iswap][m] = sendnum[iswap][m]*size_reverse;
if (m == 0) reverse_recv_offset[iswap][0] = 0;
else reverse_recv_offset[iswap][m] =
reverse_recv_offset[iswap][m-1] + sendnum[iswap][m-1];
}
rmaxswap = 0;
for (m = 0; m < nrecvproc[iswap]; m++) {
rmaxswap += recvnum[iswap][m];
size_forward_recv[iswap][m] = recvnum[iswap][m]*size_forward;
size_reverse_send[iswap][m] = recvnum[iswap][m]*size_reverse;
size_reverse_recv[iswap][m] = sendnum[iswap][m]*size_reverse;
if (m == 0) {
firstrecv[iswap][0] = atom->nlocal + atom->nghost;
forward_recv_offset[iswap][0] = 0;
@ -695,14 +952,13 @@ void CommTiled::borders()
if (ghost_velocity) {
if (sendother[iswap]) {
for (m = 0; m < nrecv; m++)
MPI_Irecv(&buf_recv[forward_recv_offset[iswap][m]],
MPI_Irecv(&buf_recv[size_border*forward_recv_offset[iswap][m]],
recvnum[iswap][m]*size_border,
MPI_DOUBLE,recvproc[iswap][m],0,world,&requests[m]);
for (m = 0; m < nsend; m++) {
n = avec->pack_border_vel(sendnum[iswap][m],sendlist[iswap][m],
buf_send,pbc_flag[iswap][m],pbc[iswap][m]);
MPI_Send(buf_send,n*size_border,MPI_DOUBLE,
sendproc[iswap][m],0,world);
MPI_Send(buf_send,n,MPI_DOUBLE,sendproc[iswap][m],0,world);
}
}
@ -718,21 +974,21 @@ void CommTiled::borders()
for (m = 0; m < nrecv; m++) {
MPI_Waitany(nrecv,requests,&irecv,&status);
avec->unpack_border(recvnum[iswap][irecv],firstrecv[iswap][irecv],
&buf_recv[forward_recv_offset[iswap][irecv]]);
&buf_recv[size_border*
forward_recv_offset[iswap][irecv]]);
}
}
} else {
if (sendother[iswap]) {
for (m = 0; m < nrecv; m++)
MPI_Irecv(&buf_recv[forward_recv_offset[iswap][m]],
MPI_Irecv(&buf_recv[size_border*forward_recv_offset[iswap][m]],
recvnum[iswap][m]*size_border,
MPI_DOUBLE,recvproc[iswap][m],0,world,&requests[m]);
for (m = 0; m < nsend; m++) {
n = avec->pack_border(sendnum[iswap][m],sendlist[iswap][m],
buf_send,pbc_flag[iswap][m],pbc[iswap][m]);
MPI_Send(buf_send,n*size_border,MPI_DOUBLE,
sendproc[iswap][m],0,world);
MPI_Send(buf_send,n,MPI_DOUBLE,sendproc[iswap][m],0,world);
}
}
@ -748,7 +1004,8 @@ void CommTiled::borders()
for (m = 0; m < nrecv; m++) {
MPI_Waitany(nrecv,requests,&irecv,&status);
avec->unpack_border(recvnum[iswap][irecv],firstrecv[iswap][irecv],
&buf_recv[forward_recv_offset[iswap][irecv]]);
&buf_recv[size_border*
forward_recv_offset[iswap][irecv]]);
}
}
}
@ -769,6 +1026,29 @@ void CommTiled::borders()
// reset global->local map
if (map_style) atom->map_set();
// DEBUG
/*
MPI_Barrier(world);
for (i = 0; i < nswap; i++) {
if (nsendproc[i] == 1)
printf("BORDERS SEND %d %d: nsend %d snum0 %d\n",
i,me,nsendproc[i],sendnum[i][0]);
else
printf("BORDERS SEND %d %d: nsend %d snums %d %d\n",
i,me,nsendproc[i],sendnum[i][0],sendnum[i][1]);
if (nrecvproc[i] == 1)
printf("BORDERS RECV %d %d: nrecv %d rnum0 %d\n",
i,me,nrecvproc[i],recvnum[i][0]);
else
printf("BORDERS RECV %d %d: nrecv %d rnums %d %d\n",
i,me,nrecvproc[i],recvnum[i][0],recvnum[i][1]);
}
MPI_Barrier(world);
*/
}
// NOTE: remaining forward/reverse methods still need to be updated
@ -1280,6 +1560,47 @@ void CommTiled::box_other_tiled(int idim, int iswap,
else hi[2] = boxhi[2];
}
/* ----------------------------------------------------------------------
return 1 if proc's box touches me, else 0
procneigh stores 6 procs that touch me
------------------------------------------------------------------------- */
int CommTiled::box_touch_brick(int proc, int idim, int iswap)
{
if (procneigh[idim][iswap] == proc) return 1;
return 0;
}
/* ----------------------------------------------------------------------
return 1 if proc's box touches me, else 0
------------------------------------------------------------------------- */
int CommTiled::box_touch_tiled(int proc, int idim, int iswap)
{
// sending to left
// only touches if proc hi = my lo, or if proc hi = boxhi and my lo = boxlo
if (iswap == 0) {
if (rcbinfo[proc].mysplit[idim][1] == rcbinfo[me].mysplit[idim][0])
return 1;
else if (rcbinfo[proc].mysplit[idim][1] == 1.0 &&
rcbinfo[me].mysplit[idim][0] == 0.0)
return 1;
// sending to right
// only touches if proc lo = my hi, or if proc lo = boxlo and my hi = boxhi
} else {
if (rcbinfo[proc].mysplit[idim][0] == rcbinfo[me].mysplit[idim][1])
return 1;
else if (rcbinfo[proc].mysplit[idim][0] == 0.0 &&
rcbinfo[me].mysplit[idim][1] == 1.0)
return 1;
}
return 0;
}
/* ----------------------------------------------------------------------
realloc the size of the send buffer as needed with BUFFACTOR and bufextra
if flag = 1, realloc
@ -1371,6 +1692,19 @@ void CommTiled::allocate_swap(int n)
grow_swap_send(i,DELTA_PROCS,0);
grow_swap_recv(i,DELTA_PROCS);
}
nesendproc = new int[n];
nerecvproc = new int[n];
neprocmax = new int[n];
esendproc = new int*[n];
erecvproc = new int*[n];
for (int i = 0; i < n; i++) {
esendproc[i] = erecvproc[i] = NULL;
neprocmax[i] = DELTA_PROCS;
esendproc[i] = new int[DELTA_PROCS];
erecvproc[i] = new int[DELTA_PROCS];
}
}
/* ----------------------------------------------------------------------
@ -1425,14 +1759,6 @@ void CommTiled::grow_swap_recv(int i, int n)
delete [] size_reverse_send[i];
size_reverse_send[i] = new int[n];
if (n > maxreqstat) {
maxreqstat = n;
delete [] requests;
delete [] statuses;
requests = new MPI_Request[n];
statuses = new MPI_Status[n];
}
}
/* ----------------------------------------------------------------------
@ -1488,6 +1814,18 @@ void CommTiled::deallocate_swap(int n)
delete [] statuses;
delete [] nprocmax;
delete [] nesendproc;
delete [] nerecvproc;
delete [] neprocmax;
for (int i = 0; i < n; i++) {
delete [] esendproc[i];
delete [] erecvproc[i];
}
delete [] esendproc;
delete [] erecvproc;
}
/* ----------------------------------------------------------------------

View File

@ -54,6 +54,9 @@ class CommTiled : public Comm {
int size_border; // # of datums in forward border comm
int nswap; // # of swaps to perform = 2*dim
// forward/reverse comm info, proc lists include self
int *nsendproc,*nrecvproc; // # of procs to send/recv to/from in each swap
int *sendother; // 1 if send to any other proc in each swap
int *sendself; // 1 if send to self in each swap
@ -74,6 +77,12 @@ class CommTiled : public Comm {
double ***sendbox; // bounding box of atoms to send per swap/proc
// exchange comm info, proc lists do not include self
int *nesendproc,*nerecvproc; // # of procs to send/recv to/from in each swap
int *neprocmax; // current max # of send procs for each swap
int **esendproc,**erecvproc; // proc to send/recv to/from per swap/proc
double *buf_send; // send buffer for all comm
double *buf_recv; // recv buffer for all comm
int maxsend,maxrecv; // current size of send/recv buffer
@ -119,6 +128,11 @@ class CommTiled : public Comm {
void box_other_brick(int, int, int, double *, double *);
void box_other_tiled(int, int, int, double *, double *);
typedef int (CommTiled::*BoxTouchPtr)(int, int, int);
BoxTouchPtr box_touch;
int box_touch_brick(int, int, int);
int box_touch_tiled(int, int, int);
void grow_send(int, int); // reallocate send buffer
void grow_recv(int); // free/allocate recv buffer
void grow_list(int, int, int); // reallocate sendlist for one swap/proc

View File

@ -266,9 +266,12 @@ void FixBalance::rebalance()
// only needed if migrate_check() says an atom moves to far
// else allow caller's comm->exchange() to do it
//NOTE: change back to migrate_check()
if (domain->triclinic) domain->x2lamda(atom->nlocal);
if (lbstyle == BISECTION) irregular->migrate_atoms(0,sendproc);
else if (irregular->migrate_check()) irregular->migrate_atoms();
//else if (irregular->migrate_check()) irregular->migrate_atoms();
else irregular->migrate_atoms();
if (domain->triclinic) domain->lamda2x(atom->nlocal);
// invoke KSpace setup_grid() to adjust to new proc sub-domains

View File

@ -158,6 +158,9 @@ void RCB::compute(int dimension, int n, double **x, double *wt,
hi[1] = bboxhi[1];
hi[2] = bboxhi[2];
cut = 0.0;
cutdim = -1;
// initialize counters
counters[0] = 0;