Files
lammps/src/npair_halffull_newtoff_trim.cpp
2022-10-24 11:08:26 -04:00

102 lines
3.0 KiB
C++

/* ----------------------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
https://www.lammps.org/, Sandia National Laboratories
LAMMPS development team: developers@lammps.org
Copyright (2003) Sandia Corporation. Under the terms of Contract
DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
certain rights in this software. This software is distributed under
the GNU General Public License.
See the README file in the top-level LAMMPS directory.
------------------------------------------------------------------------- */
#include "npair_halffull_newtoff_trim.h"
#include "atom.h"
#include "error.h"
#include "my_page.h"
#include "neigh_list.h"
using namespace LAMMPS_NS;
/* ---------------------------------------------------------------------- */
NPairHalffullNewtoffTrim::NPairHalffullNewtoffTrim(LAMMPS *lmp) : NPair(lmp) {}
/* ----------------------------------------------------------------------
build half list from full list
pair stored once if i,j are both owned and i < j
pair stored by me if j is ghost (also stored by proc owning j)
works if full list is a skip list
works for owned (non-ghost) list, also for ghost list
if ghost, also store neighbors of ghost atoms & set inum,gnum correctly
------------------------------------------------------------------------- */
void NPairHalffullNewtoffTrim::build(NeighList *list)
{
int i, j, ii, jj, n, jnum, joriginal;
int *neighptr, *jlist;
double xtmp, ytmp, ztmp;
double delx, dely, delz, rsq;
double **x = atom->x;
int *ilist = list->ilist;
int *numneigh = list->numneigh;
int **firstneigh = list->firstneigh;
MyPage<int> *ipage = list->ipage;
int *ilist_full = list->listfull->ilist;
int *numneigh_full = list->listfull->numneigh;
int **firstneigh_full = list->listfull->firstneigh;
int inum_full = list->listfull->inum;
if (list->ghost) inum_full += list->listfull->gnum;
int inum = 0;
ipage->reset();
double cutsq_custom = cutoff_custom * cutoff_custom;
// loop over atoms in full list
for (ii = 0; ii < inum_full; ii++) {
n = 0;
neighptr = ipage->vget();
// loop over parent full list
i = ilist_full[ii];
xtmp = x[i][0];
ytmp = x[i][1];
ztmp = x[i][2];
jlist = firstneigh_full[i];
jnum = numneigh_full[i];
for (jj = 0; jj < jnum; jj++) {
joriginal = jlist[jj];
j = joriginal & NEIGHMASK;
if (j > i) {
delx = xtmp - x[j][0];
dely = ytmp - x[j][1];
delz = ztmp - x[j][2];
rsq = delx * delx + dely * dely + delz * delz;
if (rsq > cutsq_custom) continue;
neighptr[n++] = joriginal;
}
}
ilist[inum++] = i;
firstneigh[i] = neighptr;
numneigh[i] = n;
ipage->vgot(n);
if (ipage->status()) error->one(FLERR, "Neighbor list overflow, boost neigh_modify one");
}
list->inum = inum;
if (list->ghost) list->gnum = list->listfull->gnum;
}