/* ---------------------------------------------------------------------- LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator https://www.lammps.org/, Sandia National Laboratories LAMMPS development team: developers@lammps.org This software is distributed under the GNU General Public License. See the README file in the top-level LAMMPS directory. ------------------------------------------------------------------------- */ /* ---------------------------------------------------------------------- Contributing author: Axel Kohlmeyer (Temple U) ------------------------------------------------------------------------- */ #include "pair_beck_omp.h" #include "atom.h" #include "comm.h" #include "force.h" #include "math_special.h" #include "neigh_list.h" #include "suffix.h" #include #include "omp_compat.h" using namespace LAMMPS_NS; using namespace MathSpecial; /* ---------------------------------------------------------------------- */ PairBeckOMP::PairBeckOMP(LAMMPS *lmp) : PairBeck(lmp), ThrOMP(lmp, THR_PAIR) { suffix_flag |= Suffix::OMP; respa_enable = 0; } /* ---------------------------------------------------------------------- */ void PairBeckOMP::compute(int eflag, int vflag) { ev_init(eflag, vflag); const int nall = atom->nlocal + atom->nghost; const int nthreads = comm->nthreads; const int inum = list->inum; #if defined(_OPENMP) #pragma omp parallel LMP_DEFAULT_NONE LMP_SHARED(eflag, vflag) #endif { int ifrom, ito, tid; loop_setup_thr(ifrom, ito, tid, inum, nthreads); ThrData *thr = fix->get_thr(tid); thr->timer(Timer::START); ev_setup_thr(eflag, vflag, nall, eatom, vatom, nullptr, thr); if (evflag) { if (eflag) { if (force->newton_pair) eval<1, 1, 1>(ifrom, ito, thr); else eval<1, 1, 0>(ifrom, ito, thr); } else { if (force->newton_pair) eval<1, 0, 1>(ifrom, ito, thr); else eval<1, 0, 0>(ifrom, ito, thr); } } else { if (force->newton_pair) eval<0, 0, 1>(ifrom, ito, thr); else eval<0, 0, 0>(ifrom, ito, thr); } thr->timer(Timer::PAIR); reduce_thr(this, eflag, vflag, thr); } // end of omp parallel region } template void PairBeckOMP::eval(int iifrom, int iito, ThrData *const thr) { int i, j, ii, jj, jnum, itype, jtype; double xtmp, ytmp, ztmp, delx, dely, delz, evdwl, fpair; double rsq, r5, force_beck, factor_lj; double r, rinv; double aaij, alphaij, betaij; double term1, term1inv, term2, term3, term4, term5, term6; int *ilist, *jlist, *numneigh, **firstneigh; evdwl = 0.0; const auto *_noalias const x = (dbl3_t *) atom->x[0]; auto *_noalias const f = (dbl3_t *) thr->get_f()[0]; int *type = atom->type; int nlocal = atom->nlocal; double *special_lj = force->special_lj; double fxtmp, fytmp, fztmp; ilist = list->ilist; numneigh = list->numneigh; firstneigh = list->firstneigh; // loop over neighbors of my atoms for (ii = iifrom; ii < iito; ++ii) { i = ilist[ii]; xtmp = x[i].x; ytmp = x[i].y; ztmp = x[i].z; itype = type[i]; jlist = firstneigh[i]; jnum = numneigh[i]; fxtmp = fytmp = fztmp = 0.0; for (jj = 0; jj < jnum; jj++) { j = jlist[jj]; factor_lj = special_lj[sbmask(j)]; j &= NEIGHMASK; delx = xtmp - x[j].x; dely = ytmp - x[j].y; delz = ztmp - x[j].z; rsq = delx * delx + dely * dely + delz * delz; jtype = type[j]; if (rsq < cutsq[itype][jtype]) { r = sqrt(rsq); r5 = rsq * rsq * r; aaij = aa[itype][jtype]; alphaij = alpha[itype][jtype]; betaij = beta[itype][jtype]; term1 = aaij * aaij + rsq; term2 = powint(term1, -5); term3 = 21.672 + 30.0 * aaij * aaij + 6.0 * rsq; term4 = alphaij + r5 * betaij; term5 = alphaij + 6.0 * r5 * betaij; rinv = 1.0 / r; force_beck = AA[itype][jtype] * exp(-1.0 * r * term4) * term5; force_beck -= BB[itype][jtype] * r * term2 * term3; fpair = factor_lj * force_beck * rinv; f[i].x += delx * fpair; f[i].y += dely * fpair; f[i].z += delz * fpair; if (NEWTON_PAIR || j < nlocal) { f[j].x -= delx * fpair; f[j].y -= dely * fpair; f[j].z -= delz * fpair; } if (EFLAG) { term6 = powint(term1, -3); term1inv = 1.0 / term1; evdwl = AA[itype][jtype] * exp(-1.0 * r * term4); evdwl -= BB[itype][jtype] * term6 * (1.0 + (2.709 + 3.0 * aaij * aaij) * term1inv); evdwl *= factor_lj; } if (EVFLAG) ev_tally_thr(this, i, j, nlocal, NEWTON_PAIR, evdwl, 0.0, fpair, delx, dely, delz, thr); } } f[i].x += fxtmp; f[i].y += fytmp; f[i].z += fztmp; } } /* ---------------------------------------------------------------------- */ double PairBeckOMP::memory_usage() { double bytes = memory_usage_thr(); bytes += PairBeck::memory_usage(); return bytes; }