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

This commit is contained in:
sjplimp 2010-04-23 18:45:31 +00:00
parent f994e8a1e9
commit d1bcc71748
4 changed files with 68 additions and 30 deletions

View File

@ -231,9 +231,10 @@ void Neighbor::granular_nsq_newton(NeighList *list)
if ((itag+jtag) % 2 == 1) continue;
} else {
if (x[j][2] < ztmp) continue;
else if (x[j][2] == ztmp && x[j][1] < ytmp) continue;
else if (x[j][2] == ztmp && x[j][1] == ytmp && x[j][0] < xtmp)
continue;
if (x[j][2] == ztmp) {
if (x[j][1] < ytmp) continue;
if (x[j][1] == ytmp && x[j][0] < xtmp) continue;
}
}
}
@ -477,11 +478,14 @@ void Neighbor::granular_bin_newton(NeighList *list)
for (j = bins[i]; j >= 0; j = bins[j]) {
if (j >= nlocal) {
if (x[j][2] < ztmp) continue;
if (x[j][2] == ztmp && x[j][1] < ytmp) continue;
if (x[j][2] == ztmp && x[j][1] == ytmp && x[j][0] < xtmp) continue;
if (exclude && exclusion(i,j,type[i],type[j],mask,molecule)) continue;
if (x[j][2] == ztmp) {
if (x[j][1] < ytmp) continue;
if (x[j][1] == ytmp && x[j][0] < xtmp) continue;
}
}
if (exclude && exclusion(i,j,type[i],type[j],mask,molecule)) continue;
delx = xtmp - x[j][0];
dely = ytmp - x[j][1];
delz = ztmp - x[j][2];
@ -579,15 +583,22 @@ void Neighbor::granular_bin_newton_tri(NeighList *list)
// loop over all atoms in bins in stencil
// pairs for atoms j "below" i are excluded
// below = lower z or (equal z and lower y) or (equal zy and <= x)
// this excludes self-self interaction
// below = lower z or (equal z and lower y) or (equal zy and lower x)
// (equal zyx and j <= i)
// latter excludes self-self interaction but allows superposed atoms
ibin = coord2bin(x[i]);
for (k = 0; k < nstencil; k++) {
for (j = binhead[ibin+stencil[k]]; j >= 0; j = bins[j]) {
if (x[j][2] < ztmp) continue;
if (x[j][2] == ztmp && x[j][1] < ytmp) continue;
if (x[j][2] == ztmp && x[j][1] == ytmp && x[j][0] <= xtmp) continue;
if (x[j][2] == ztmp) {
if (x[j][1] < ytmp) continue;
if (x[j][1] == ytmp) {
if (x[j][0] < xtmp) continue;
if (x[j][0] == xtmp && j <= i) continue;
}
}
if (exclude && exclusion(i,j,type[i],type[j],mask,molecule)) continue;
delx = xtmp - x[j][0];

View File

@ -173,8 +173,10 @@ void Neighbor::half_bin_newton(NeighList *list)
for (j = bins[i]; j >= 0; j = bins[j]) {
if (j >= nlocal) {
if (x[j][2] < ztmp) continue;
if (x[j][2] == ztmp && x[j][1] < ytmp) continue;
if (x[j][2] == ztmp && x[j][1] == ytmp && x[j][0] < xtmp) continue;
if (x[j][2] == ztmp) {
if (x[j][1] < ytmp) continue;
if (x[j][1] == ytmp && x[j][0] < xtmp) continue;
}
}
jtype = type[j];
@ -282,15 +284,21 @@ void Neighbor::half_bin_newton_tri(NeighList *list)
// loop over all atoms in bins in stencil
// pairs for atoms j "below" i are excluded
// below = lower z or (equal z and lower y) or (equal zy and <= x)
// this excludes self-self interaction
// below = lower z or (equal z and lower y) or (equal zy and lower x)
// (equal zyx and j <= i)
// latter excludes self-self interaction but allows superposed atoms
ibin = coord2bin(x[i]);
for (k = 0; k < nstencil; k++) {
for (j = binhead[ibin+stencil[k]]; j >= 0; j = bins[j]) {
if (x[j][2] < ztmp) continue;
if (x[j][2] == ztmp && x[j][1] < ytmp) continue;
if (x[j][2] == ztmp && x[j][1] == ytmp && x[j][0] <= xtmp) continue;
if (x[j][2] == ztmp) {
if (x[j][1] < ytmp) continue;
if (x[j][1] == ytmp) {
if (x[j][0] < xtmp) continue;
if (x[j][0] == xtmp && j <= i) continue;
}
}
jtype = type[j];
if (exclude && exclusion(i,j,itype,jtype,mask,molecule)) continue;

View File

@ -184,8 +184,10 @@ void Neighbor::half_multi_newton(NeighList *list)
for (j = bins[i]; j >= 0; j = bins[j]) {
if (j >= nlocal) {
if (x[j][2] < ztmp) continue;
if (x[j][2] == ztmp && x[j][1] < ytmp) continue;
if (x[j][2] == ztmp && x[j][1] == ytmp && x[j][0] < xtmp) continue;
if (x[j][2] == ztmp) {
if (x[j][1] < ytmp) continue;
if (x[j][1] == ytmp && x[j][0] < xtmp) continue;
}
}
jtype = type[j];
@ -304,7 +306,10 @@ void Neighbor::half_multi_newton_tri(NeighList *list)
// loop over all atoms in bins, including self, in stencil
// skip if i,j neighbor cutoff is less than bin distance
// bins below self are excluded from stencil
// pairs for atoms j below i are excluded
// pairs for atoms j "below" i are excluded
// below = lower z or (equal z and lower y) or (equal zy and lower x)
// (equal zyx and j <= i)
// latter excludes self-self interaction but allows superposed atoms
ibin = coord2bin(x[i]);
s = stencil_multi[itype];
@ -316,8 +321,13 @@ void Neighbor::half_multi_newton_tri(NeighList *list)
jtype = type[j];
if (cutsq[jtype] < distsq[k]) continue;
if (x[j][2] < ztmp) continue;
if (x[j][2] == ztmp && x[j][1] < ytmp) continue;
if (x[j][2] == ztmp && x[j][1] == ytmp && x[j][0] <= xtmp) continue;
if (x[j][2] == ztmp) {
if (x[j][1] < ytmp) continue;
if (x[j][1] == ytmp) {
if (x[j][0] < xtmp) continue;
if (x[j][0] == xtmp && j <= i) continue;
}
}
if (exclude && exclusion(i,j,itype,jtype,mask,molecule)) continue;

View File

@ -276,9 +276,10 @@ void Neighbor::respa_nsq_newton(NeighList *list)
if ((itag+jtag) % 2 == 1) continue;
} else {
if (x[j][2] < ztmp) continue;
else if (x[j][2] == ztmp && x[j][1] < ytmp) continue;
else if (x[j][2] == ztmp && x[j][1] == ytmp && x[j][0] < xtmp)
continue;
else if (x[j][2] == ztmp) {
if (x[j][1] < ytmp) continue;
if (x[j][1] == ytmp && x[j][0] < xtmp) continue;
}
}
}
@ -607,8 +608,10 @@ void Neighbor::respa_bin_newton(NeighList *list)
for (j = bins[i]; j >= 0; j = bins[j]) {
if (j >= nlocal) {
if (x[j][2] < ztmp) continue;
if (x[j][2] == ztmp && x[j][1] < ytmp) continue;
if (x[j][2] == ztmp && x[j][1] == ytmp && x[j][0] < xtmp) continue;
if (x[j][2] == ztmp) {
if (x[j][1] < ytmp) continue;
if (x[j][1] == ytmp && x[j][0] < xtmp) continue;
}
}
jtype = type[j];
@ -799,15 +802,21 @@ void Neighbor::respa_bin_newton_tri(NeighList *list)
// loop over all atoms in bins in stencil
// pairs for atoms j "below" i are excluded
// below = lower z or (equal z and lower y) or (equal zy and <= x)
// this excludes self-self interaction
// below = lower z or (equal z and lower y) or (equal zy and lower x)
// (equal zyx and j <= i)
// latter excludes self-self interaction but allows superposed atoms
ibin = coord2bin(x[i]);
for (k = 0; k < nstencil; k++) {
for (j = binhead[ibin+stencil[k]]; j >= 0; j = bins[j]) {
if (x[j][2] < ztmp) continue;
if (x[j][2] == ztmp && x[j][1] < ytmp) continue;
if (x[j][2] == ztmp && x[j][1] == ytmp && x[j][0] <= xtmp) continue;
if (x[j][2] == ztmp) {
if (x[j][1] < ytmp) continue;
if (x[j][1] == ytmp) {
if (x[j][0] < xtmp) continue;
if (x[j][0] == xtmp && j <= i) continue;
}
}
jtype = type[j];
if (exclude && exclusion(i,j,itype,jtype,mask,molecule)) continue;