forked from lijiext/lammps
394 lines
9.8 KiB
Plaintext
394 lines
9.8 KiB
Plaintext
/* ----------------------------------------------------------------------
|
|
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: Inderaj Bains (NVIDIA), ibains@nvidia.com
|
|
------------------------------------------------------------------------- */
|
|
|
|
#ifndef LJE_GPU_KERNEL
|
|
#define LJE_GPU_KERNEL
|
|
|
|
#ifdef NV_KERNEL
|
|
|
|
#include "nv_kernel_def.h"
|
|
texture<float4> pos_tex;
|
|
|
|
#ifdef _DOUBLE_DOUBLE
|
|
__inline double4 fetch_pos(const int& i, const double4 *pos)
|
|
{
|
|
return pos[i];
|
|
}
|
|
#else
|
|
__inline float4 fetch_pos(const int& i, const float4 *pos)
|
|
{
|
|
return tex1Dfetch(pos_tex, i);
|
|
}
|
|
#endif
|
|
|
|
#else
|
|
|
|
#pragma OPENCL EXTENSION cl_khr_fp64: enable
|
|
#define GLOBAL_ID_X get_global_id(0)
|
|
#define THREAD_ID_X get_local_id(0)
|
|
#define BLOCK_ID_X get_group_id(0)
|
|
#define BLOCK_SIZE_X get_local_size(0)
|
|
#define __syncthreads() barrier(CLK_LOCAL_MEM_FENCE)
|
|
#define __inline inline
|
|
|
|
#define fetch_pos(i,y) x_[i]
|
|
#define BLOCK_PAIR 64
|
|
#define MAX_SHARED_TYPES 8
|
|
|
|
#endif
|
|
|
|
#ifdef _DOUBLE_DOUBLE
|
|
#define numtyp double
|
|
#define numtyp2 double2
|
|
#define numtyp4 double4
|
|
#define acctyp double
|
|
#define acctyp4 double4
|
|
#endif
|
|
|
|
#ifdef _SINGLE_DOUBLE
|
|
#define numtyp float
|
|
#define numtyp2 float2
|
|
#define numtyp4 float4
|
|
#define acctyp double
|
|
#define acctyp4 double4
|
|
#endif
|
|
|
|
#ifndef numtyp
|
|
#define numtyp float
|
|
#define numtyp2 float2
|
|
#define numtyp4 float4
|
|
#define acctyp float
|
|
#define acctyp4 float4
|
|
#endif
|
|
|
|
#define SBBITS 30
|
|
#define NEIGHMASK 0x3FFFFFFF
|
|
__inline int sbmask(int j) { return j >> SBBITS & 3; }
|
|
|
|
__kernel void kernel_pair(__global numtyp4 *x_, __global numtyp4 *lj1,
|
|
__global numtyp4* lj3, const int lj_types,
|
|
__global numtyp *sp_lj_in, __global int *dev_nbor,
|
|
__global int *dev_packed, __global acctyp4 *ans,
|
|
__global acctyp *engv, const int eflag,
|
|
const int vflag, const int inum, const int nall,
|
|
const int nbor_pitch, const int t_per_atom) {
|
|
int tid=THREAD_ID_X;
|
|
int ii=mul24((int)BLOCK_ID_X,(int)(BLOCK_SIZE_X)/t_per_atom);
|
|
ii+=tid/t_per_atom;
|
|
int offset=tid%t_per_atom;
|
|
|
|
__local numtyp sp_lj[4];
|
|
sp_lj[0]=sp_lj_in[0];
|
|
sp_lj[1]=sp_lj_in[1];
|
|
sp_lj[2]=sp_lj_in[2];
|
|
sp_lj[3]=sp_lj_in[3];
|
|
|
|
acctyp energy=(acctyp)0;
|
|
acctyp4 f;
|
|
f.x=(acctyp)0;
|
|
f.y=(acctyp)0;
|
|
f.z=(acctyp)0;
|
|
acctyp virial[6];
|
|
for (int i=0; i<6; i++)
|
|
virial[i]=(acctyp)0;
|
|
|
|
if (ii<inum) {
|
|
__global int *nbor=dev_nbor+ii;
|
|
int i=*nbor;
|
|
nbor+=nbor_pitch;
|
|
int numj=*nbor;
|
|
nbor+=nbor_pitch;
|
|
|
|
int n_stride;
|
|
__global int *list_end;
|
|
if (dev_nbor==dev_packed) {
|
|
list_end=nbor+mul24(numj,nbor_pitch);
|
|
nbor+=mul24(offset,nbor_pitch);
|
|
n_stride=mul24(t_per_atom,nbor_pitch);
|
|
} else {
|
|
nbor=dev_packed+*nbor;
|
|
list_end=nbor+numj;
|
|
n_stride=t_per_atom;
|
|
nbor+=offset;
|
|
}
|
|
|
|
numtyp4 ix=fetch_pos(i,x_); //x_[i];
|
|
int itype=ix.w;
|
|
|
|
numtyp factor_lj;
|
|
for ( ; nbor<list_end; nbor+=n_stride) {
|
|
|
|
int j=*nbor;
|
|
factor_lj = sp_lj[sbmask(j)];
|
|
j &= NEIGHMASK;
|
|
|
|
numtyp4 jx=fetch_pos(j,x_); //x_[j];
|
|
int jtype=jx.w;
|
|
|
|
// Compute r12
|
|
numtyp delx = ix.x-jx.x;
|
|
numtyp dely = ix.y-jx.y;
|
|
numtyp delz = ix.z-jx.z;
|
|
numtyp r2inv = delx*delx+dely*dely+delz*delz;
|
|
|
|
int mtype=itype*lj_types+jtype;
|
|
if (r2inv<lj1[mtype].z) {
|
|
numtyp r = sqrt(r2inv);
|
|
numtyp rshift = r - lj1[mtype].w;
|
|
numtyp rshiftsq = rshift*rshift;
|
|
r2inv = (numtyp) 1.0/rshiftsq;
|
|
numtyp r6inv = r2inv*r2inv*r2inv;
|
|
numtyp force = r6inv*(lj1[mtype].x*r6inv-lj1[mtype].y);
|
|
force*=factor_lj/rshift/r;
|
|
|
|
f.x+=delx*force;
|
|
f.y+=dely*force;
|
|
f.z+=delz*force;
|
|
|
|
if (eflag>0) {
|
|
numtyp e=r6inv*(lj3[mtype].x*r6inv-lj3[mtype].y);
|
|
energy+=factor_lj*(e-lj3[mtype].z);
|
|
}
|
|
if (vflag>0) {
|
|
virial[0] += delx*delx*force;
|
|
virial[1] += dely*dely*force;
|
|
virial[2] += delz*delz*force;
|
|
virial[3] += delx*dely*force;
|
|
virial[4] += delx*delz*force;
|
|
virial[5] += dely*delz*force;
|
|
}
|
|
}
|
|
|
|
} // for nbor
|
|
} // if ii
|
|
|
|
// Reduce answers
|
|
if (t_per_atom>1) {
|
|
__local acctyp red_acc[6][BLOCK_PAIR];
|
|
|
|
red_acc[0][tid]=f.x;
|
|
red_acc[1][tid]=f.y;
|
|
red_acc[2][tid]=f.z;
|
|
red_acc[3][tid]=energy;
|
|
|
|
for (unsigned int s=t_per_atom/2; s>0; s>>=1) {
|
|
if (offset < s) {
|
|
for (int r=0; r<4; r++)
|
|
red_acc[r][tid] += red_acc[r][tid+s];
|
|
}
|
|
}
|
|
|
|
f.x=red_acc[0][tid];
|
|
f.y=red_acc[1][tid];
|
|
f.z=red_acc[2][tid];
|
|
energy=red_acc[3][tid];
|
|
|
|
if (vflag>0) {
|
|
for (int r=0; r<6; r++)
|
|
red_acc[r][tid]=virial[r];
|
|
|
|
for (unsigned int s=t_per_atom/2; s>0; s>>=1) {
|
|
if (offset < s) {
|
|
for (int r=0; r<6; r++)
|
|
red_acc[r][tid] += red_acc[r][tid+s];
|
|
}
|
|
}
|
|
|
|
for (int r=0; r<6; r++)
|
|
virial[r]=red_acc[r][tid];
|
|
}
|
|
}
|
|
|
|
// Store answers
|
|
if (ii<inum && offset==0) {
|
|
__global acctyp *ap1=engv+ii;
|
|
if (eflag>0) {
|
|
*ap1=energy;
|
|
ap1+=inum;
|
|
}
|
|
if (vflag>0) {
|
|
for (int i=0; i<6; i++) {
|
|
*ap1=virial[i];
|
|
ap1+=inum;
|
|
}
|
|
}
|
|
ans[ii]=f;
|
|
} // if ii
|
|
}
|
|
|
|
__kernel void kernel_pair_fast(__global numtyp4 *x_, __global numtyp4 *lj1_in,
|
|
__global numtyp4* lj3_in,
|
|
__global numtyp* sp_lj_in,
|
|
__global int *dev_nbor, __global int *dev_packed,
|
|
__global acctyp4 *ans, __global acctyp *engv,
|
|
const int eflag, const int vflag, const int inum,
|
|
const int nall, const int nbor_pitch,
|
|
const int t_per_atom) {
|
|
int tid=THREAD_ID_X;
|
|
int ii=mul24((int)BLOCK_ID_X,(int)(BLOCK_SIZE_X)/t_per_atom);
|
|
ii+=tid/t_per_atom;
|
|
int offset=tid%t_per_atom;
|
|
|
|
__local numtyp4 lj1[MAX_SHARED_TYPES*MAX_SHARED_TYPES];
|
|
__local numtyp4 lj3[MAX_SHARED_TYPES*MAX_SHARED_TYPES];
|
|
__local numtyp sp_lj[4];
|
|
if (tid<4)
|
|
sp_lj[tid]=sp_lj_in[tid];
|
|
if (tid<MAX_SHARED_TYPES*MAX_SHARED_TYPES) {
|
|
lj1[tid]=lj1_in[tid];
|
|
if (eflag>0)
|
|
lj3[tid]=lj3_in[tid];
|
|
}
|
|
|
|
acctyp energy=(acctyp)0;
|
|
acctyp4 f;
|
|
f.x=(acctyp)0;
|
|
f.y=(acctyp)0;
|
|
f.z=(acctyp)0;
|
|
acctyp virial[6];
|
|
for (int i=0; i<6; i++)
|
|
virial[i]=(numtyp)0;
|
|
|
|
__syncthreads();
|
|
|
|
if (ii<inum) {
|
|
__global int *nbor=dev_nbor+ii;
|
|
int i=*nbor;
|
|
nbor+=nbor_pitch;
|
|
int numj=*nbor;
|
|
nbor+=nbor_pitch;
|
|
|
|
int n_stride;
|
|
__global int *list_end;
|
|
if (dev_nbor==dev_packed) {
|
|
list_end=nbor+mul24(numj,nbor_pitch);
|
|
nbor+=mul24(offset,nbor_pitch);
|
|
n_stride=mul24(t_per_atom,nbor_pitch);
|
|
} else {
|
|
nbor=dev_packed+*nbor;
|
|
list_end=nbor+numj;
|
|
n_stride=t_per_atom;
|
|
nbor+=offset;
|
|
}
|
|
|
|
numtyp4 ix=fetch_pos(i,x_); //x_[i];
|
|
int iw=ix.w;
|
|
int itype=mul24((int)MAX_SHARED_TYPES,iw);
|
|
|
|
numtyp factor_lj;
|
|
for ( ; nbor<list_end; nbor+=n_stride) {
|
|
|
|
int j=*nbor;
|
|
factor_lj = sp_lj[sbmask(j)];
|
|
j &= NEIGHMASK;
|
|
|
|
numtyp4 jx=fetch_pos(j,x_); //x_[j];
|
|
int mtype=itype+jx.w;
|
|
|
|
// Compute r12
|
|
numtyp delx = ix.x-jx.x;
|
|
numtyp dely = ix.y-jx.y;
|
|
numtyp delz = ix.z-jx.z;
|
|
numtyp r2inv = delx*delx+dely*dely+delz*delz;
|
|
|
|
if (r2inv<lj1[mtype].z) {
|
|
numtyp r = sqrt(r2inv);
|
|
numtyp rshift = r - lj1[mtype].w;
|
|
numtyp rshiftsq = rshift*rshift;
|
|
r2inv = 1.0/rshiftsq;
|
|
numtyp r6inv = r2inv*r2inv*r2inv;
|
|
numtyp force = r6inv*(lj1[mtype].x*r6inv-lj1[mtype].y);
|
|
force*=factor_lj/rshift/r;
|
|
|
|
f.x+=delx*force;
|
|
f.y+=dely*force;
|
|
f.z+=delz*force;
|
|
|
|
if (eflag>0) {
|
|
numtyp e=r6inv*(lj3[mtype].x*r6inv-lj3[mtype].y);
|
|
energy+=factor_lj*(e-lj3[mtype].z);
|
|
}
|
|
if (vflag>0) {
|
|
virial[0] += delx*delx*force;
|
|
virial[1] += dely*dely*force;
|
|
virial[2] += delz*delz*force;
|
|
virial[3] += delx*dely*force;
|
|
virial[4] += delx*delz*force;
|
|
virial[5] += dely*delz*force;
|
|
}
|
|
}
|
|
|
|
} // for nbor
|
|
} // if ii
|
|
|
|
// Reduce answers
|
|
if (t_per_atom>1) {
|
|
__local acctyp red_acc[6][BLOCK_PAIR];
|
|
|
|
red_acc[0][tid]=f.x;
|
|
red_acc[1][tid]=f.y;
|
|
red_acc[2][tid]=f.z;
|
|
red_acc[3][tid]=energy;
|
|
|
|
for (unsigned int s=t_per_atom/2; s>0; s>>=1) {
|
|
if (offset < s) {
|
|
for (int r=0; r<4; r++)
|
|
red_acc[r][tid] += red_acc[r][tid+s];
|
|
}
|
|
}
|
|
|
|
f.x=red_acc[0][tid];
|
|
f.y=red_acc[1][tid];
|
|
f.z=red_acc[2][tid];
|
|
energy=red_acc[3][tid];
|
|
|
|
if (vflag>0) {
|
|
for (int r=0; r<6; r++)
|
|
red_acc[r][tid]=virial[r];
|
|
|
|
for (unsigned int s=t_per_atom/2; s>0; s>>=1) {
|
|
if (offset < s) {
|
|
for (int r=0; r<6; r++)
|
|
red_acc[r][tid] += red_acc[r][tid+s];
|
|
}
|
|
}
|
|
|
|
for (int r=0; r<6; r++)
|
|
virial[r]=red_acc[r][tid];
|
|
}
|
|
}
|
|
|
|
// Store answers
|
|
if (ii<inum && offset==0) {
|
|
__global acctyp *ap1=engv+ii;
|
|
if (eflag>0) {
|
|
*ap1=energy;
|
|
ap1+=inum;
|
|
}
|
|
if (vflag>0) {
|
|
for (int i=0; i<6; i++) {
|
|
*ap1=virial[i];
|
|
ap1+=inum;
|
|
}
|
|
}
|
|
ans[ii]=f;
|
|
} // if ii*/
|
|
}
|
|
|
|
#endif
|
|
|