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_skip_kokkos.h"
16 
17 #include "atom_kokkos.h"
18 #include "atom_masks.h"
19 #include "atom_vec.h"
20 #include "neigh_list_kokkos.h"
21 
22 using namespace LAMMPS_NS;
23 
24 /* ---------------------------------------------------------------------- */
25 
26 template<class DeviceType>
NPairSkipKokkos(LAMMPS * lmp)27 NPairSkipKokkos<DeviceType>::NPairSkipKokkos(LAMMPS *lmp) : NPair(lmp) {
28   atomKK = (AtomKokkos *) atom;
29   execution_space = ExecutionSpaceFromDevice<DeviceType>::space;
30   d_inum = typename AT::t_int_scalar("npair_skip:inum");
31 }
32 
33 /* ----------------------------------------------------------------------
34    build skip list for subset of types from parent list
35    works for half and full lists
36    works for owned (non-ghost) list, also for ghost list
37    iskip and ijskip flag which atom types and type pairs to skip
38    if ghost, also store neighbors of ghost atoms & set inum,gnum correctly
39 ------------------------------------------------------------------------- */
40 
41 template<class DeviceType>
build(NeighList * list)42 void NPairSkipKokkos<DeviceType>::build(NeighList *list)
43 {
44   atomKK->sync(execution_space,TYPE_MASK);
45   type = atomKK->k_type.view<DeviceType>();
46   nlocal = atom->nlocal;
47 
48 
49   NeighListKokkos<DeviceType>* k_list_skip = static_cast<NeighListKokkos<DeviceType>*>(list->listskip);
50   d_ilist_skip = k_list_skip->d_ilist;
51   d_numneigh_skip = k_list_skip->d_numneigh;
52   d_neighbors_skip = k_list_skip->d_neighbors;
53 
54   num_skip = list->listskip->inum;
55   if (list->ghost) num_skip += list->listskip->gnum;
56 
57   NeighListKokkos<DeviceType>* k_list = static_cast<NeighListKokkos<DeviceType>*>(list);
58   k_list->maxneighs = k_list_skip->maxneighs; // simple, but could be made more memory efficient
59   k_list->grow(atom->nmax);
60   d_ilist = k_list->d_ilist;
61   d_numneigh = k_list->d_numneigh;
62   d_neighbors = k_list->d_neighbors;
63 
64   int ntypes = atom->ntypes;
65 
66   k_iskip = DAT::tdual_int_1d("npair_skip:iskip",ntypes+1);
67   k_ijskip = DAT::tdual_int_2d("npair_skip:ijskip",ntypes+1,ntypes+1);
68   d_iskip = k_iskip.view<DeviceType>();
69   d_ijskip = k_ijskip.view<DeviceType>();
70 
71   for (int itype = 1; itype <= ntypes; itype++) {
72     k_iskip.h_view(itype) = list->iskip[itype];
73     for (int jtype = 1; jtype <= ntypes; jtype++) {
74       k_ijskip.h_view(itype,jtype) = list->ijskip[itype][jtype];
75     }
76   }
77   k_iskip.modify<LMPHostType>();
78   k_ijskip.modify<LMPHostType>();
79 
80   k_iskip.sync<DeviceType>();
81   k_ijskip.sync<DeviceType>();
82 
83   // loop over atoms in other list
84   // skip I atom entirely if iskip is set for type[I]
85   // skip I,J pair if ijskip is set for type[I],type[J]
86 
87   copymode = 1;
88   Kokkos::parallel_scan(Kokkos::RangePolicy<DeviceType, TagNPairSkipCompute>(0,num_skip),*this);
89 
90   auto h_inum = Kokkos::create_mirror_view(d_inum);
91   Kokkos::deep_copy(h_inum,d_inum);
92   const int inum = h_inum();
93   list->inum = inum;
94   if (list->ghost) {
95     int num = 0;
96     Kokkos::parallel_reduce(Kokkos::RangePolicy<DeviceType, TagNPairSkipCountLocal>(0,inum),*this,num);
97     list->inum = num;
98     list->gnum = inum - num;
99   }
100   copymode = 0;
101 }
102 
103 template<class DeviceType>
104 KOKKOS_INLINE_FUNCTION
operator ()(TagNPairSkipCompute,const int & ii,int & inum,const bool & final) const105 void NPairSkipKokkos<DeviceType>::operator()(TagNPairSkipCompute, const int &ii, int &inum, const bool &final) const {
106 
107   const int i = d_ilist_skip(ii);
108   const int itype = type(i);
109 
110   if (!d_iskip(itype)) {
111 
112     if (final) {
113 
114       int n = 0;
115 
116       // loop over parent non-skip list
117 
118       const int jnum = d_numneigh_skip(i);
119       const AtomNeighbors neighbors_i = AtomNeighbors(&d_neighbors(i,0),d_numneigh(i),
120                                                     &d_neighbors(i,1)-&d_neighbors(i,0));
121 
122       for (int jj = 0; jj < jnum; jj++) {
123         const int joriginal = d_neighbors_skip(i,jj);
124         int j = joriginal & NEIGHMASK;
125         if (d_ijskip(itype,type(j))) continue;
126         neighbors_i(n++) = joriginal;
127       }
128 
129       d_numneigh(i) = n;
130       d_ilist(inum) = i;
131     }
132 
133     inum++;
134   }
135 
136   if (final) {
137     if (ii == num_skip-1)
138       d_inum() = inum;
139   }
140 }
141 
142 template<class DeviceType>
143 KOKKOS_INLINE_FUNCTION
operator ()(TagNPairSkipCountLocal,const int & i,int & num) const144 void NPairSkipKokkos<DeviceType>::operator()(TagNPairSkipCountLocal, const int &i, int &num) const {
145   if (d_ilist[i] < nlocal) num++;
146 }
147 
148 
149 namespace LAMMPS_NS {
150 template class NPairSkipKokkos<LMPDeviceType>;
151 #ifdef LMP_KOKKOS_GPU
152 template class NPairSkipKokkos<LMPHostType>;
153 #endif
154 }
155