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_respa_nsq_newton.h"
16 #include "neigh_list.h"
17 #include "atom.h"
18 #include "atom_vec.h"
19 #include "group.h"
20 #include "molecule.h"
21 #include "domain.h"
22 #include "my_page.h"
23 #include "error.h"
24 
25 using namespace LAMMPS_NS;
26 
27 /* ---------------------------------------------------------------------- */
28 
NPairHalfRespaNsqNewton(LAMMPS * lmp)29 NPairHalfRespaNsqNewton::NPairHalfRespaNsqNewton(LAMMPS *lmp) : NPair(lmp) {}
30 
31 /* ----------------------------------------------------------------------
32    multiple respa lists
33    N^2 / 2 search for neighbor pairs with full Newton's 3rd law
34    pair added to list if atoms i and j are both owned and i < j
35    if j is ghost only me or other proc adds pair
36    decision based on itag,jtag tests
37 ------------------------------------------------------------------------- */
38 
build(NeighList * list)39 void NPairHalfRespaNsqNewton::build(NeighList *list)
40 {
41   int i,j,n,itype,jtype,itag,jtag,n_inner,n_middle,bitmask;
42   int imol,iatom,moltemplate;
43   tagint tagprev;
44   double xtmp,ytmp,ztmp,delx,dely,delz,rsq;
45   int *neighptr,*neighptr_inner,*neighptr_middle;
46 
47   double **x = atom->x;
48   int *type = atom->type;
49   int *mask = atom->mask;
50   tagint *tag = atom->tag;
51   tagint *molecule = atom->molecule;
52   tagint **special = atom->special;
53   int **nspecial = atom->nspecial;
54   int nlocal = atom->nlocal;
55   int nall = nlocal + atom->nghost;
56   if (includegroup) {
57     nlocal = atom->nfirst;
58     bitmask = group->bitmask[includegroup];
59   }
60 
61   int *molindex = atom->molindex;
62   int *molatom = atom->molatom;
63   Molecule **onemols = atom->avec->onemols;
64   if (molecular == Atom::TEMPLATE) moltemplate = 1;
65   else moltemplate = 0;
66 
67   int *ilist = list->ilist;
68   int *numneigh = list->numneigh;
69   int **firstneigh = list->firstneigh;
70   MyPage<int> *ipage = list->ipage;
71 
72   int *ilist_inner = list->ilist_inner;
73   int *numneigh_inner = list->numneigh_inner;
74   int **firstneigh_inner = list->firstneigh_inner;
75   MyPage<int> *ipage_inner = list->ipage_inner;
76 
77   int *ilist_middle,*numneigh_middle,**firstneigh_middle;
78   MyPage<int> *ipage_middle;
79   int respamiddle = list->respamiddle;
80   if (respamiddle) {
81     ilist_middle = list->ilist_middle;
82     numneigh_middle = list->numneigh_middle;
83     firstneigh_middle = list->firstneigh_middle;
84     ipage_middle = list->ipage_middle;
85   }
86 
87   int inum = 0;
88   int which = 0;
89   int minchange = 0;
90   ipage->reset();
91   ipage_inner->reset();
92   if (respamiddle) ipage_middle->reset();
93 
94   for (i = 0; i < nlocal; i++) {
95     n = n_inner = 0;
96     neighptr = ipage->vget();
97     neighptr_inner = ipage_inner->vget();
98     if (respamiddle) {
99       n_middle = 0;
100       neighptr_middle = ipage_middle->vget();
101     }
102 
103     itag = tag[i];
104     itype = type[i];
105     xtmp = x[i][0];
106     ytmp = x[i][1];
107     ztmp = x[i][2];
108     if (moltemplate) {
109       imol = molindex[i];
110       iatom = molatom[i];
111       tagprev = tag[i] - iatom - 1;
112     }
113 
114     // loop over remaining atoms, owned and ghost
115 
116     for (j = i+1; j < nall; j++) {
117       if (includegroup && !(mask[j] & bitmask)) continue;
118 
119       if (j >= nlocal) {
120         jtag = tag[j];
121         if (itag > jtag) {
122           if ((itag+jtag) % 2 == 0) continue;
123         } else if (itag < jtag) {
124           if ((itag+jtag) % 2 == 1) continue;
125         } else {
126           if (x[j][2] < ztmp) continue;
127           if (x[j][2] == ztmp) {
128             if (x[j][1] < ytmp) continue;
129             if (x[j][1] == ytmp && x[j][0] < xtmp) continue;
130           }
131         }
132       }
133 
134       jtype = type[j];
135       if (exclude && exclusion(i,j,itype,jtype,mask,molecule)) continue;
136 
137       delx = xtmp - x[j][0];
138       dely = ytmp - x[j][1];
139       delz = ztmp - x[j][2];
140       rsq = delx*delx + dely*dely + delz*delz;
141 
142       if (rsq <= cutneighsq[itype][jtype]) {
143         if (molecular != Atom::ATOMIC) {
144           if (!moltemplate)
145             which = find_special(special[i],nspecial[i],tag[j]);
146           else if (imol >= 0)
147             which = find_special(onemols[imol]->special[iatom],
148                                  onemols[imol]->nspecial[iatom],
149                                  tag[j]-tagprev);
150           else which = 0;
151           if (which == 0) neighptr[n++] = j;
152           else if ((minchange = domain->minimum_image_check(delx,dely,delz)))
153             neighptr[n++] = j;
154           else if (which > 0) neighptr[n++] = j ^ (which << SBBITS);
155         } else neighptr[n++] = j;
156 
157         if (rsq < cut_inner_sq) {
158           if (which == 0) neighptr_inner[n_inner++] = j;
159           else if (minchange) neighptr_inner[n_inner++] = j;
160           else if (which > 0) neighptr_inner[n_inner++] = j ^ (which << SBBITS);
161         }
162 
163         if (respamiddle &&
164             rsq < cut_middle_sq && rsq > cut_middle_inside_sq) {
165           if (which == 0) neighptr_middle[n_middle++] = j;
166           else if (minchange) neighptr_middle[n_middle++] = j;
167           else if (which > 0)
168             neighptr_middle[n_middle++] = j ^ (which << SBBITS);
169         }
170       }
171     }
172 
173     ilist[inum] = i;
174     firstneigh[i] = neighptr;
175     numneigh[i] = n;
176     ipage->vgot(n);
177     if (ipage->status())
178       error->one(FLERR,"Neighbor list overflow, boost neigh_modify one");
179 
180     ilist_inner[inum] = i;
181     firstneigh_inner[i] = neighptr_inner;
182     numneigh_inner[i] = n_inner;
183     ipage_inner->vgot(n_inner);
184     if (ipage_inner->status())
185       error->one(FLERR,"Neighbor list overflow, boost neigh_modify one");
186 
187     if (respamiddle) {
188       ilist_middle[inum] = i;
189       firstneigh_middle[i] = neighptr_middle;
190       numneigh_middle[i] = n_middle;
191       ipage_middle->vgot(n_middle);
192       if (ipage_middle->status())
193         error->one(FLERR,"Neighbor list overflow, boost neigh_modify one");
194     }
195 
196     inum++;
197   }
198 
199   list->inum = inum;
200   list->inum_inner = inum;
201   if (respamiddle) list->inum_middle = inum;
202 }
203