499 lines
16 KiB
C++
499 lines
16 KiB
C++
/* -*- c++ -*- ----------------------------------------------------------
|
|
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
|
|
http://lammps.sandia.gov, Sandia National Laboratories
|
|
Steve Plimpton, sjplimp@sandia.gov
|
|
|
|
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.
|
|
------------------------------------------------------------------------- */
|
|
|
|
/* ----------------------------------------------------------------------
|
|
Contributing author: Rodrigo Canales (RWTH Aachen University)
|
|
------------------------------------------------------------------------- */
|
|
|
|
#include <math.h>
|
|
#include "pair_buck_intel.h"
|
|
#include "atom.h"
|
|
#include "comm.h"
|
|
#include "force.h"
|
|
#include "group.h"
|
|
#include "neighbor.h"
|
|
#include "neigh_list.h"
|
|
#include "neigh_request.h"
|
|
#include "math_const.h"
|
|
#include "memory.h"
|
|
#include "suffix.h"
|
|
#include "force.h"
|
|
#include "modify.h"
|
|
|
|
using namespace LAMMPS_NS;
|
|
using namespace MathConst;
|
|
|
|
#define C_FORCE_T typename ForceConst<flt_t>::c_force_t
|
|
#define C_ENERGY_T typename ForceConst<flt_t>::c_energy_t
|
|
|
|
PairBuckIntel::PairBuckIntel(LAMMPS *lmp) : PairBuck(lmp)
|
|
{
|
|
suffix_flag |= Suffix::INTEL;
|
|
}
|
|
|
|
PairBuckIntel::~PairBuckIntel()
|
|
{
|
|
}
|
|
|
|
void PairBuckIntel::compute(int eflag, int vflag)
|
|
{
|
|
if (fix->precision()==FixIntel::PREC_MODE_MIXED)
|
|
compute<float,double>(eflag, vflag, fix->get_mixed_buffers(),
|
|
force_const_single);
|
|
else if (fix->precision()==FixIntel::PREC_MODE_DOUBLE)
|
|
compute<double,double>(eflag, vflag, fix->get_double_buffers(),
|
|
force_const_double);
|
|
else
|
|
compute<float,float>(eflag, vflag, fix->get_single_buffers(),
|
|
force_const_single);
|
|
|
|
fix->balance_stamp();
|
|
vflag_fdotr = 0;
|
|
}
|
|
|
|
template <class flt_t, class acc_t>
|
|
void PairBuckIntel::compute(int eflag, int vflag,
|
|
IntelBuffers<flt_t,acc_t> *buffers,
|
|
const ForceConst<flt_t> &fc)
|
|
{
|
|
if (eflag || vflag) {
|
|
ev_setup(eflag,vflag);
|
|
} else evflag = vflag_fdotr = 0;
|
|
|
|
const int inum = list->inum;
|
|
const int nthreads = comm->nthreads;
|
|
const int host_start = fix->host_start_pair();
|
|
const int offload_end = fix->offload_end_pair();
|
|
const int ago = neighbor->ago;
|
|
|
|
if (ago != 0 && fix->separate_buffers() == 0) {
|
|
fix->start_watch(TIME_PACK);
|
|
#if defined(_OPENMP)
|
|
#pragma omp parallel default(none) shared(eflag,vflag,buffers,fc)
|
|
#endif
|
|
{
|
|
int ifrom, ito, tid;
|
|
IP_PRE_omp_range_id_align(ifrom, ito, tid, atom->nlocal + atom->nghost,
|
|
nthreads, sizeof(ATOM_T));
|
|
buffers->thr_pack(ifrom,ito,ago);
|
|
}
|
|
fix->stop_watch(TIME_PACK);
|
|
}
|
|
|
|
if (evflag || vflag_fdotr) {
|
|
int ovflag = 0;
|
|
if (vflag_fdotr) ovflag = 2;
|
|
else if (vflag) ovflag = 1;
|
|
if (eflag) {
|
|
if (force->newton_pair) {
|
|
eval<1,1,1>(1, ovflag, buffers, fc, 0, offload_end);
|
|
eval<1,1,1>(0, ovflag, buffers, fc, host_start, inum);
|
|
} else {
|
|
eval<1,1,0>(1, ovflag, buffers, fc, 0, offload_end);
|
|
eval<1,1,0>(0, ovflag, buffers, fc, host_start, inum);
|
|
}
|
|
} else {
|
|
if (force->newton_pair) {
|
|
eval<1,0,1>(1, ovflag, buffers, fc, 0, offload_end);
|
|
eval<1,0,1>(0, ovflag, buffers, fc, host_start, inum);
|
|
} else {
|
|
eval<1,0,0>(1, ovflag, buffers, fc, 0, offload_end);
|
|
eval<1,0,0>(0, ovflag, buffers, fc, host_start, inum);
|
|
}
|
|
}
|
|
} else {
|
|
if (force->newton_pair) {
|
|
eval<0,0,1>(1, 0, buffers, fc, 0, offload_end);
|
|
eval<0,0,1>(0, 0, buffers, fc, host_start, inum);
|
|
} else {
|
|
eval<0,0,0>(1, 0, buffers, fc, 0, offload_end);
|
|
eval<0,0,0>(0, 0, buffers, fc, host_start, inum);
|
|
}
|
|
}
|
|
}
|
|
|
|
/* ---------------------------------------------------------------------- */
|
|
|
|
template <int EVFLAG, int EFLAG, int NEWTON_PAIR, class flt_t, class acc_t>
|
|
void PairBuckIntel::eval(const int offload, const int vflag,
|
|
IntelBuffers<flt_t,acc_t> *buffers,
|
|
const ForceConst<flt_t> &fc,
|
|
const int astart, const int aend)
|
|
{
|
|
const int inum = aend - astart;
|
|
if (inum == 0) return;
|
|
int nlocal, nall, minlocal;
|
|
fix->get_buffern(offload, nlocal, nall, minlocal);
|
|
const int ago = neighbor->ago;
|
|
IP_PRE_pack_separate_buffers(fix, buffers, ago, offload, nlocal, nall);
|
|
|
|
ATOM_T * _noalias const x = buffers->get_x(offload);
|
|
|
|
const int * _noalias const numneigh = list->numneigh;
|
|
const int * _noalias const cnumneigh = buffers->cnumneigh(list);
|
|
const int * _noalias const firstneigh = buffers->firstneigh(list);
|
|
|
|
const flt_t * _noalias const special_lj = fc.special_lj;
|
|
const C_FORCE_T * _noalias const c_force = fc.c_force[0];
|
|
const C_ENERGY_T * _noalias const c_energy = fc.c_energy[0];
|
|
|
|
const int ntypes = atom->ntypes + 1;
|
|
const int eatom = this->eflag_atom;
|
|
|
|
// Determine how much data to transfer
|
|
int x_size, q_size, f_stride, ev_size, separate_flag;
|
|
IP_PRE_get_transfern(ago, NEWTON_PAIR, EVFLAG, EFLAG, vflag,
|
|
buffers, offload, fix, separate_flag,
|
|
x_size, q_size, ev_size, f_stride);
|
|
|
|
int tc;
|
|
FORCE_T * _noalias f_start;
|
|
acc_t * _noalias ev_global;
|
|
IP_PRE_get_buffers(offload, buffers, fix, tc, f_start, ev_global);
|
|
|
|
const int nthreads = tc;
|
|
#ifdef _LMP_INTEL_OFFLOAD
|
|
int *overflow = fix->get_off_overflow_flag();
|
|
double *timer_compute = fix->off_watch_pair();
|
|
// Redeclare as local variables for offload
|
|
|
|
if (offload) fix->start_watch(TIME_OFFLOAD_LATENCY);
|
|
#pragma offload target(mic:_cop) if(offload) \
|
|
in(special_lj:length(0) alloc_if(0) free_if(0)) \
|
|
in(c_force, c_energy:length(0) alloc_if(0) free_if(0)) \
|
|
in(firstneigh:length(0) alloc_if(0) free_if(0)) \
|
|
in(cnumneigh:length(0) alloc_if(0) free_if(0)) \
|
|
in(numneigh:length(0) alloc_if(0) free_if(0)) \
|
|
in(x:length(x_size) alloc_if(0) free_if(0)) \
|
|
in(overflow:length(0) alloc_if(0) free_if(0)) \
|
|
in(astart,nthreads,inum,nall,ntypes,vflag,eatom) \
|
|
in(f_stride,nlocal,minlocal,separate_flag,offload) \
|
|
out(f_start:length(f_stride) alloc_if(0) free_if(0)) \
|
|
out(ev_global:length(ev_size) alloc_if(0) free_if(0)) \
|
|
out(timer_compute:length(1) alloc_if(0) free_if(0)) \
|
|
signal(f_start)
|
|
#endif
|
|
{
|
|
#if defined(__MIC__) && defined(_LMP_INTEL_OFFLOAD)
|
|
*timer_compute = MIC_Wtime();
|
|
#endif
|
|
|
|
IP_PRE_repack_for_offload(NEWTON_PAIR, separate_flag, nlocal, nall,
|
|
f_stride, x, 0);
|
|
|
|
acc_t oevdwl, ov0, ov1, ov2, ov3, ov4, ov5;
|
|
if (EVFLAG) {
|
|
oevdwl = (acc_t)0;
|
|
if (vflag) ov0 = ov1 = ov2 = ov3 = ov4 = ov5 = (acc_t)0;
|
|
}
|
|
|
|
// loop over neighbors of my atoms
|
|
#if defined(_OPENMP)
|
|
#pragma omp parallel default(none) \
|
|
shared(f_start,f_stride,nlocal,nall,minlocal) \
|
|
reduction(+:oevdwl,ov0,ov1,ov2,ov3,ov4,ov5)
|
|
#endif
|
|
{
|
|
int iifrom, iito, tid;
|
|
IP_PRE_omp_range_id(iifrom, iito, tid, inum, nthreads);
|
|
iifrom += astart;
|
|
iito += astart;
|
|
|
|
FORCE_T * _noalias const f = f_start - minlocal + (tid * f_stride);
|
|
memset(f + minlocal, 0, f_stride * sizeof(FORCE_T));
|
|
|
|
for (int i = iifrom; i < iito; ++i) {
|
|
const int itype = x[i].w;
|
|
|
|
const int ptr_off = itype * ntypes;
|
|
const C_FORCE_T * _noalias const c_forcei = c_force + ptr_off;
|
|
const C_ENERGY_T * _noalias const c_energyi = c_energy + ptr_off;
|
|
const int * _noalias const jlist = firstneigh + cnumneigh[i];
|
|
const int jnum = numneigh[i];
|
|
|
|
acc_t fxtmp,fytmp,fztmp,fwtmp;
|
|
acc_t sevdwl, sv0, sv1, sv2, sv3, sv4, sv5;
|
|
|
|
const flt_t xtmp = x[i].x;
|
|
const flt_t ytmp = x[i].y;
|
|
const flt_t ztmp = x[i].z;
|
|
fxtmp = fytmp = fztmp = (acc_t)0;
|
|
if (EVFLAG) {
|
|
if (EFLAG) fwtmp = sevdwl = (acc_t)0;
|
|
if (vflag==1) sv0 = sv1 = sv2 = sv3 = sv4 = sv5 = (acc_t)0;
|
|
}
|
|
|
|
#if defined(LMP_SIMD_COMPILER)
|
|
#pragma vector aligned
|
|
#pragma simd reduction(+:fxtmp, fytmp, fztmp, fwtmp, sevdwl, \
|
|
sv0, sv1, sv2, sv3, sv4, sv5)
|
|
#endif
|
|
for (int jj = 0; jj < jnum; jj++) {
|
|
|
|
flt_t forcebuck, evdwl;
|
|
forcebuck = evdwl = (flt_t)0.0;
|
|
|
|
const int sbindex = jlist[jj] >> SBBITS & 3;
|
|
const int j = jlist[jj] & NEIGHMASK;
|
|
|
|
const flt_t delx = xtmp - x[j].x;
|
|
const flt_t dely = ytmp - x[j].y;
|
|
const flt_t delz = ztmp - x[j].z;
|
|
const int jtype = x[j].w;
|
|
const flt_t rsq = delx * delx + dely * dely + delz * delz;
|
|
const flt_t r = sqrt(rsq);
|
|
const flt_t r2inv = (flt_t)1.0 / rsq;
|
|
|
|
#ifdef INTEL_VMASK
|
|
if (rsq < c_forcei[jtype].cutsq) {
|
|
#endif
|
|
const flt_t r6inv = r2inv * r2inv * r2inv;
|
|
const flt_t rexp = exp(-r * c_forcei[jtype].rhoinv);
|
|
forcebuck = r * rexp * c_forcei[jtype].buck1 -
|
|
r6inv * c_forcei[jtype].buck2;
|
|
|
|
#ifndef INTEL_VMASK
|
|
if (rsq > c_forcei[jtype].cutsq)
|
|
forcebuck =(flt_t)0.0;
|
|
#endif
|
|
if (EFLAG) {
|
|
evdwl = rexp * c_energyi[jtype].a -
|
|
r6inv * c_energyi[jtype].c -
|
|
c_energyi[jtype].offset;
|
|
|
|
#ifndef INTEL_VMASK
|
|
if (rsq > c_forcei[jtype].cutsq)
|
|
evdwl =(flt_t)0.0;
|
|
#endif
|
|
}
|
|
|
|
if (sbindex) {
|
|
const flt_t factor_lj = special_lj[sbindex];
|
|
forcebuck *= factor_lj;
|
|
if (EFLAG)
|
|
evdwl *= factor_lj;
|
|
}
|
|
const flt_t fpair = forcebuck * r2inv;
|
|
fxtmp += delx * fpair;
|
|
fytmp += dely * fpair;
|
|
fztmp += delz * fpair;
|
|
if (NEWTON_PAIR || j < nlocal) {
|
|
f[j].x -= delx * fpair;
|
|
f[j].y -= dely * fpair;
|
|
f[j].z -= delz * fpair;
|
|
}
|
|
|
|
if (EVFLAG) {
|
|
flt_t ev_pre = (flt_t)0;
|
|
if (NEWTON_PAIR || i < nlocal)
|
|
ev_pre += (flt_t)0.5;
|
|
if (NEWTON_PAIR || j < nlocal)
|
|
ev_pre += (flt_t)0.5;
|
|
|
|
if (EFLAG) {
|
|
sevdwl += ev_pre * evdwl;
|
|
if (eatom) {
|
|
if (NEWTON_PAIR || i < nlocal)
|
|
fwtmp += (flt_t)0.5 * evdwl;
|
|
if (NEWTON_PAIR || j < nlocal)
|
|
f[j].w += (flt_t)0.5 * evdwl;
|
|
}
|
|
}
|
|
IP_PRE_ev_tally_nbor(vflag, ev_pre, fpair, delx, dely, delz);
|
|
}
|
|
#ifdef INTEL_VMASK
|
|
}
|
|
#endif
|
|
} // for jj
|
|
|
|
f[i].x += fxtmp;
|
|
f[i].y += fytmp;
|
|
f[i].z += fztmp;
|
|
IP_PRE_ev_tally_atom(EVFLAG, EFLAG, vflag, f, fwtmp);
|
|
} // for ii
|
|
|
|
#ifndef _LMP_INTEL_OFFLOAD
|
|
if (vflag == 2)
|
|
#endif
|
|
{
|
|
#if defined(_OPENMP)
|
|
#pragma omp barrier
|
|
#endif
|
|
IP_PRE_fdotr_acc_force(NEWTON_PAIR, EVFLAG, EFLAG, vflag, eatom, nall,
|
|
nlocal, minlocal, nthreads, f_start, f_stride,
|
|
x, offload);
|
|
}
|
|
} // end of omp parallel region
|
|
if (EVFLAG) {
|
|
if (EFLAG) {
|
|
ev_global[0] = oevdwl;
|
|
ev_global[1] = (acc_t)0;
|
|
}
|
|
if (vflag) {
|
|
ev_global[2] = ov0;
|
|
ev_global[3] = ov1;
|
|
ev_global[4] = ov2;
|
|
ev_global[5] = ov3;
|
|
ev_global[6] = ov4;
|
|
ev_global[7] = ov5;
|
|
}
|
|
}
|
|
#if defined(__MIC__) && defined(_LMP_INTEL_OFFLOAD)
|
|
*timer_compute = MIC_Wtime() - *timer_compute;
|
|
#endif
|
|
} // end of offload region
|
|
|
|
if (offload)
|
|
fix->stop_watch(TIME_OFFLOAD_LATENCY);
|
|
else
|
|
fix->stop_watch(TIME_HOST_PAIR);
|
|
|
|
if (EVFLAG)
|
|
fix->add_result_array(f_start, ev_global, offload, eatom, 0, vflag);
|
|
else
|
|
fix->add_result_array(f_start, 0, offload);
|
|
}
|
|
|
|
void PairBuckIntel::init_style()
|
|
{
|
|
PairBuck::init_style();
|
|
neighbor->requests[neighbor->nrequest-1]->intel = 1;
|
|
|
|
int ifix = modify->find_fix("package_intel");
|
|
if (ifix < 0)
|
|
error->all(FLERR,
|
|
"The 'package intel' command is required for /intel styles");
|
|
fix = static_cast<FixIntel *>(modify->fix[ifix]);
|
|
|
|
fix->pair_init_check();
|
|
#ifdef _LMP_INTEL_OFFLOAD
|
|
_cop = fix->coprocessor_number();
|
|
#endif
|
|
|
|
if (fix->precision() == FixIntel::PREC_MODE_MIXED)
|
|
pack_force_const(force_const_single, fix->get_mixed_buffers());
|
|
else if (fix->precision() == FixIntel::PREC_MODE_DOUBLE)
|
|
pack_force_const(force_const_double, fix->get_double_buffers());
|
|
else
|
|
pack_force_const(force_const_single, fix->get_single_buffers());
|
|
}
|
|
|
|
template <class flt_t, class acc_t>
|
|
void PairBuckIntel::pack_force_const(ForceConst<flt_t> &fc,
|
|
IntelBuffers<flt_t,acc_t> *buffers)
|
|
{
|
|
int tp1 = atom->ntypes + 1;
|
|
|
|
fc.set_ntypes(tp1, memory, _cop);
|
|
buffers->set_ntypes(tp1);
|
|
flt_t **cutneighsq = buffers->get_cutneighsq();
|
|
|
|
// Repeat cutsq calculation because done after call to init_style
|
|
double cut, cutneigh;
|
|
for (int i = 1; i <= atom->ntypes; i++) {
|
|
for (int j = i; j <= atom->ntypes; j++) {
|
|
if (setflag[i][j] != 0 || (setflag[i][i] != 0 && setflag[j][j] != 0)) {
|
|
cut = init_one(i, j);
|
|
cutneigh = cut + neighbor->skin;
|
|
cutsq[i][j] = cutsq[j][i] = cut*cut;
|
|
cutneighsq[i][j] = cutneighsq[j][i] = cutneigh * cutneigh;
|
|
}
|
|
}
|
|
}
|
|
|
|
for (int i = 0; i < 4; i++) {
|
|
fc.special_lj[i] = force->special_lj[i];
|
|
fc.special_lj[0] = 1.0;
|
|
}
|
|
|
|
for (int i = 0; i < tp1; i++) {
|
|
for (int j = 0; j < tp1; j++) {
|
|
fc.c_force[i][j].buck1 = buck1[i][j];
|
|
fc.c_force[i][j].buck2 = buck2[i][j];
|
|
fc.c_force[i][j].rhoinv = rhoinv[i][j];
|
|
fc.c_force[i][j].cutsq = cutsq[i][j];
|
|
fc.c_energy[i][j].a = a[i][j];
|
|
fc.c_energy[i][j].c = c[i][j];
|
|
fc.c_energy[i][j].offset = offset[i][j];
|
|
}
|
|
}
|
|
|
|
#ifdef _LMP_INTEL_OFFLOAD
|
|
if (_cop < 0) return;
|
|
flt_t * special_lj = fc.special_lj;
|
|
C_FORCE_T * c_force = fc.c_force[0];
|
|
C_ENERGY_T * c_energy = fc.c_energy[0];
|
|
flt_t * ocutneighsq = cutneighsq[0];
|
|
int tp1sq = tp1 * tp1;
|
|
#pragma offload_transfer target(mic:_cop) \
|
|
in(special_lj: length(4) alloc_if(0) free_if(0)) \
|
|
in(c_force, c_energy: length(tp1sq) alloc_if(0) free_if(0)) \
|
|
in(ocutneighsq: length(tp1sq) alloc_if(0) free_if(0))
|
|
#endif
|
|
}
|
|
|
|
/* ---------------------------------------------------------------------- */
|
|
|
|
template <class flt_t>
|
|
void PairBuckIntel::ForceConst<flt_t>::set_ntypes(const int ntypes,
|
|
Memory *memory,
|
|
const int cop) {
|
|
if ( (ntypes != _ntypes ) ) {
|
|
if (_ntypes > 0) {
|
|
#ifdef _LMP_INTEL_OFFLOAD
|
|
flt_t * ospecial_lj = special_lj;
|
|
c_force_t * oc_force = c_force[0];
|
|
c_energy_t * oc_energy = c_energy[0];
|
|
|
|
if (ospecial_lj != NULL && oc_force != NULL &&
|
|
oc_energy != NULL &&
|
|
_cop >= 0) {
|
|
#pragma offload_transfer target(mic:cop) \
|
|
nocopy(ospecial_lj: alloc_if(0) free_if(1)) \
|
|
nocopy(oc_force, oc_energy: alloc_if(0) free_if(1))
|
|
}
|
|
#endif
|
|
|
|
_memory->destroy(c_force);
|
|
_memory->destroy(c_energy);
|
|
|
|
}
|
|
if (ntypes > 0) {
|
|
_cop = cop;
|
|
memory->create(c_force,ntypes,ntypes,"fc.c_force");
|
|
memory->create(c_energy,ntypes,ntypes,"fc.c_energy");
|
|
|
|
|
|
#ifdef _LMP_INTEL_OFFLOAD
|
|
flt_t * ospecial_lj = special_lj;
|
|
c_force_t * oc_force = c_force[0];
|
|
c_energy_t * oc_energy = c_energy[0];
|
|
int tp1sq = ntypes*ntypes;
|
|
if (ospecial_lj != NULL && oc_force != NULL &&
|
|
oc_energy != NULL &&
|
|
cop >= 0) {
|
|
#pragma offload_transfer target(mic:cop) \
|
|
nocopy(ospecial_lj: length(4) alloc_if(1) free_if(0)) \
|
|
nocopy(oc_force: length(tp1sq) alloc_if(1) free_if(0)) \
|
|
nocopy(oc_energy: length(tp1sq) alloc_if(1) free_if(0))
|
|
|
|
}
|
|
#endif
|
|
}
|
|
}
|
|
_ntypes=ntypes;
|
|
_memory=memory;
|
|
}
|
|
|
|
|