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_full_bin_atomonly.h"
16
17 #include "atom.h"
18 #include "error.h"
19 #include "my_page.h"
20 #include "neigh_list.h"
21
22 using namespace LAMMPS_NS;
23
24 /* ---------------------------------------------------------------------- */
25
NPairFullBinAtomonly(LAMMPS * lmp)26 NPairFullBinAtomonly::NPairFullBinAtomonly(LAMMPS *lmp) : NPair(lmp) {}
27
28 /* ----------------------------------------------------------------------
29 binned neighbor list construction for all neighbors
30 every neighbor pair appears in list of both atoms i and j
31 ------------------------------------------------------------------------- */
32
build(NeighList * list)33 void NPairFullBinAtomonly::build(NeighList *list)
34 {
35 int i,j,k,n,itype,jtype,ibin;
36 double xtmp,ytmp,ztmp,delx,dely,delz,rsq;
37 int *neighptr;
38
39 double **x = atom->x;
40 int *type = atom->type;
41 int *mask = atom->mask;
42 tagint *molecule = atom->molecule;
43 int nlocal = atom->nlocal;
44 if (includegroup) nlocal = atom->nfirst;
45
46 int *ilist = list->ilist;
47 int *numneigh = list->numneigh;
48 int **firstneigh = list->firstneigh;
49 MyPage<int> *ipage = list->ipage;
50
51 int inum = 0;
52 ipage->reset();
53
54 for (i = 0; i < nlocal; i++) {
55 n = 0;
56 neighptr = ipage->vget();
57
58 itype = type[i];
59 xtmp = x[i][0];
60 ytmp = x[i][1];
61 ztmp = x[i][2];
62
63 // loop over all atoms in surrounding bins in stencil including self
64 // skip i = j
65
66 ibin = atom2bin[i];
67
68 for (k = 0; k < nstencil; k++) {
69 for (j = binhead[ibin+stencil[k]]; j >= 0; j = bins[j]) {
70 if (i == j) continue;
71
72 jtype = type[j];
73 if (exclude && exclusion(i,j,itype,jtype,mask,molecule)) continue;
74
75 delx = xtmp - x[j][0];
76 dely = ytmp - x[j][1];
77 delz = ztmp - x[j][2];
78 rsq = delx*delx + dely*dely + delz*delz;
79
80 if (rsq <= cutneighsq[itype][jtype]) neighptr[n++] = j;
81 }
82 }
83
84 ilist[inum++] = i;
85 firstneigh[i] = neighptr;
86 numneigh[i] = n;
87 ipage->vgot(n);
88 if (ipage->status())
89 error->one(FLERR,"Neighbor list overflow, boost neigh_modify one");
90 }
91
92 list->inum = inum;
93 list->gnum = 0;
94 }
95