1 // clang-format off
2 /* ----------------------------------------------------------------------
3    LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
4    https://www.lammps.org/, Sandia National Laboratories
5    Steve Plimpton, sjplimp@sandia.gov
6 
7    Copyright (2003) Sandia Corporation.  Under the terms of Contract
8    DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
9    certain rights in this software.  This software is distributed under
10    the GNU General Public License.
11 
12    See the README file in the top-level LAMMPS directory.
13 ------------------------------------------------------------------------- */
14 
15 #include "npair_half_bin_newton.h"
16 #include "neigh_list.h"
17 #include "atom.h"
18 #include "atom_vec.h"
19 #include "molecule.h"
20 #include "domain.h"
21 #include "my_page.h"
22 #include "error.h"
23 
24 using namespace LAMMPS_NS;
25 
26 /* ---------------------------------------------------------------------- */
27 
NPairHalfBinNewton(LAMMPS * lmp)28 NPairHalfBinNewton::NPairHalfBinNewton(LAMMPS *lmp) : NPair(lmp) {}
29 
30 /* ----------------------------------------------------------------------
31    binned neighbor list construction with full Newton's 3rd law
32    each owned atom i checks its own bin and other bins in Newton stencil
33    every pair stored exactly once by some processor
34 ------------------------------------------------------------------------- */
35 
build(NeighList * list)36 void NPairHalfBinNewton::build(NeighList *list)
37 {
38   int i,j,k,n,itype,jtype,ibin,which,imol,iatom,moltemplate;
39   tagint tagprev;
40   double xtmp,ytmp,ztmp,delx,dely,delz,rsq;
41   int *neighptr;
42 
43   double **x = atom->x;
44   int *type = atom->type;
45   int *mask = atom->mask;
46   tagint *tag = atom->tag;
47   tagint *molecule = atom->molecule;
48   tagint **special = atom->special;
49   int **nspecial = atom->nspecial;
50   int nlocal = atom->nlocal;
51   if (includegroup) nlocal = atom->nfirst;
52 
53   int *molindex = atom->molindex;
54   int *molatom = atom->molatom;
55   Molecule **onemols = atom->avec->onemols;
56   if (molecular == Atom::TEMPLATE) moltemplate = 1;
57   else moltemplate = 0;
58 
59   int *ilist = list->ilist;
60   int *numneigh = list->numneigh;
61   int **firstneigh = list->firstneigh;
62   MyPage<int> *ipage = list->ipage;
63 
64   int inum = 0;
65   ipage->reset();
66 
67   for (i = 0; i < nlocal; i++) {
68     n = 0;
69     neighptr = ipage->vget();
70 
71     itype = type[i];
72     xtmp = x[i][0];
73     ytmp = x[i][1];
74     ztmp = x[i][2];
75     if (moltemplate) {
76       imol = molindex[i];
77       iatom = molatom[i];
78       tagprev = tag[i] - iatom - 1;
79     }
80 
81     // loop over rest of atoms in i's bin, ghosts are at end of linked list
82     // if j is owned atom, store it, since j is beyond i in linked list
83     // if j is ghost, only store if j coords are "above and to the right" of i
84 
85     for (j = bins[i]; j >= 0; j = bins[j]) {
86       if (j >= nlocal) {
87         if (x[j][2] < ztmp) continue;
88         if (x[j][2] == ztmp) {
89           if (x[j][1] < ytmp) continue;
90           if (x[j][1] == ytmp && x[j][0] < xtmp) continue;
91         }
92       }
93 
94       jtype = type[j];
95       if (exclude && exclusion(i,j,itype,jtype,mask,molecule)) continue;
96 
97       delx = xtmp - x[j][0];
98       dely = ytmp - x[j][1];
99       delz = ztmp - x[j][2];
100       rsq = delx*delx + dely*dely + delz*delz;
101 
102       if (rsq <= cutneighsq[itype][jtype]) {
103         if (molecular != Atom::ATOMIC) {
104           if (!moltemplate)
105             which = find_special(special[i],nspecial[i],tag[j]);
106           else if (imol >= 0)
107             which = find_special(onemols[imol]->special[iatom],
108                                  onemols[imol]->nspecial[iatom],
109                                  tag[j]-tagprev);
110           else which = 0;
111           if (which == 0) neighptr[n++] = j;
112           else if (domain->minimum_image_check(delx,dely,delz))
113             neighptr[n++] = j;
114           else if (which > 0) neighptr[n++] = j ^ (which << SBBITS);
115           // OLD: if (which >= 0) neighptr[n++] = j ^ (which << SBBITS);
116         } else neighptr[n++] = j;
117       }
118     }
119 
120     // loop over all atoms in other bins in stencil, store every pair
121 
122     ibin = atom2bin[i];
123     for (k = 0; k < nstencil; k++) {
124       for (j = binhead[ibin+stencil[k]]; j >= 0; j = bins[j]) {
125         jtype = type[j];
126         if (exclude && exclusion(i,j,itype,jtype,mask,molecule)) continue;
127 
128         delx = xtmp - x[j][0];
129         dely = ytmp - x[j][1];
130         delz = ztmp - x[j][2];
131         rsq = delx*delx + dely*dely + delz*delz;
132 
133         if (rsq <= cutneighsq[itype][jtype]) {
134           if (molecular != Atom::ATOMIC) {
135             if (!moltemplate)
136               which = find_special(special[i],nspecial[i],tag[j]);
137             else if (imol >= 0)
138               which = find_special(onemols[imol]->special[iatom],
139                                    onemols[imol]->nspecial[iatom],
140                                    tag[j]-tagprev);
141             else which = 0;
142             if (which == 0) neighptr[n++] = j;
143             else if (domain->minimum_image_check(delx,dely,delz))
144               neighptr[n++] = j;
145             else if (which > 0) neighptr[n++] = j ^ (which << SBBITS);
146             // OLD: if (which >= 0) neighptr[n++] = j ^ (which << SBBITS);
147           } else neighptr[n++] = j;
148         }
149       }
150     }
151 
152     ilist[inum++] = i;
153     firstneigh[i] = neighptr;
154     numneigh[i] = n;
155     ipage->vgot(n);
156     if (ipage->status())
157       error->one(FLERR,"Neighbor list overflow, boost neigh_modify one");
158   }
159 
160   list->inum = inum;
161 }
162