/* ---------------------------------------------------------------------- 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 authors: Mike Brown (ORNL), brownw@ornl.gov ------------------------------------------------------------------------- */ #ifndef LJC_GPU_KERNEL #define LJC_GPU_KERNEL #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 #ifdef NV_KERNEL #include "nv_kernel_def.h" texture pos_tex; texture q_tex; #ifdef _DOUBLE_DOUBLE __inline double4 fetch_pos(const int& i, const double4 *pos) { return pos[i]; } __inline double fetch_q(const int& i, const double *q) { return q[i]; } #else __inline float4 fetch_pos(const int& i, const float4 *pos) { return tex1Dfetch(pos_tex, i); } __inline float fetch_q(const int& i, const float *q) { return tex1Dfetch(q_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 fetch_q(i,y) q_[i] #define BLOCK_PAIR 64 #define MAX_SHARED_TYPES 8 #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, __global numtyp *q_ , __global numtyp *cutsq, const numtyp qqrd2e, 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[8]; 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]; sp_lj[4]=sp_lj_in[4]; sp_lj[5]=sp_lj_in[5]; sp_lj[6]=sp_lj_in[6]; sp_lj[7]=sp_lj_in[7]; acctyp energy=(acctyp)0; acctyp e_coul=(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 (ii0) { e_coul += forcecoul; if (rsq < lj1[mtype].z) { 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; red_acc[4][tid]=e_coul; for (unsigned int s=t_per_atom/2; s>0; s>>=1) { if (offset < s) { for (int r=0; r<5; 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]; e_coul=red_acc[4][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 (ii0) { *ap1=energy; ap1+=inum; *ap1=e_coul; 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, __global numtyp *q_ , __global numtyp *_cutsq, const numtyp qqrd2e, 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 cutsq[MAX_SHARED_TYPES*MAX_SHARED_TYPES]; __local numtyp sp_lj[8]; if (tid<8) sp_lj[tid]=sp_lj_in[tid]; if (tid0) lj3[tid]=lj3_in[tid]; } acctyp energy=(acctyp)0; acctyp e_coul=(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; __syncthreads(); if (ii0) { e_coul += forcecoul; if (rsq < lj1[mtype].z) { 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; red_acc[4][tid]=e_coul; for (unsigned int s=t_per_atom/2; s>0; s>>=1) { if (offset < s) { for (int r=0; r<5; 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]; e_coul=red_acc[4][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 (ii0) { *ap1=energy; ap1+=inum; *ap1=e_coul; ap1+=inum; } if (vflag>0) { for (int i=0; i<6; i++) { *ap1=virial[i]; ap1+=inum; } } ans[ii]=f; } // if ii*/ } #endif