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