Implement GPU pair style lj/cut/tip4p/long/gpu

Source code, Makefiles and Install for GPU-accelerated TIP4P pair style.
It is implemented as a part of the standard GPU package.
The style is compatible with the standard  lj/cut/tip4p/long.
Also, this commit modifies "atom.h" just to
add a getter for variable 'max_same'.
This commit is contained in:
Vsevak 2019-11-10 02:38:58 +03:00
parent 90729ebe25
commit 64bdc59623
10 changed files with 1287 additions and 4 deletions

View File

@ -82,7 +82,8 @@ OBJS = $(OBJ_DIR)/lal_atom.o $(OBJ_DIR)/lal_ans.o \
$(OBJ_DIR)/lal_lj_expand_coul_long.o $(OBJ_DIR)/lal_lj_expand_coul_long_ext.o \
$(OBJ_DIR)/lal_coul_long_cs.o $(OBJ_DIR)/lal_coul_long_cs_ext.o \
$(OBJ_DIR)/lal_born_coul_long_cs.o $(OBJ_DIR)/lal_born_coul_long_cs_ext.o \
$(OBJ_DIR)/lal_born_coul_wolf_cs.o $(OBJ_DIR)/lal_born_coul_wolf_cs_ext.o
$(OBJ_DIR)/lal_born_coul_wolf_cs.o $(OBJ_DIR)/lal_born_coul_wolf_cs_ext.o \
$(OBJ_DIR)/lal_lj_tip4p_long.o $(OBJ_DIR)/lal_lj_tip4p_long_ext.o
CBNS = $(OBJ_DIR)/device.cubin $(OBJ_DIR)/device_cubin.h \
$(OBJ_DIR)/atom.cubin $(OBJ_DIR)/atom_cubin.h \
@ -143,7 +144,8 @@ CBNS = $(OBJ_DIR)/device.cubin $(OBJ_DIR)/device_cubin.h \
$(OBJ_DIR)/lj_expand_coul_long.cubin $(OBJ_DIR)/lj_expand_coul_long_cubin.h \
$(OBJ_DIR)/coul_long_cs.cubin $(OBJ_DIR)/coul_long_cs_cubin.h \
$(OBJ_DIR)/born_coul_long_cs.cubin $(OBJ_DIR)/born_coul_long_cs_cubin.h \
$(OBJ_DIR)/born_coul_wolf_cs.cubin $(OBJ_DIR)/born_coul_wolf_cs_cubin.h
$(OBJ_DIR)/born_coul_wolf_cs.cubin $(OBJ_DIR)/born_coul_wolf_cs_cubin.h \
$(OBJ_DIR)/lj_tip4p_long.cubin $(OBJ_DIR)/lj_tip4p_long_cubin.h
all: $(OBJ_DIR) $(GPU_LIB) $(EXECS)
@ -297,6 +299,18 @@ $(OBJ_DIR)/lal_lj.o: $(ALL_H) lal_lj.h lal_lj.cpp $(OBJ_DIR)/lj_cubin.h $(OBJ_DI
$(OBJ_DIR)/lal_lj_ext.o: $(ALL_H) lal_lj.h lal_lj_ext.cpp lal_base_atomic.h
$(CUDR) -o $@ -c lal_lj_ext.cpp -I$(OBJ_DIR)
$(OBJ_DIR)/lj_tip4p_long.cubin: lal_lj_tip4p_long.cu lal_precision.h lal_preprocessor.h
$(CUDA) --cubin -DNV_KERNEL -o $@ lal_lj_tip4p_long.cu
$(OBJ_DIR)/lj_tip4p_long_cubin.h: $(OBJ_DIR)/lj_tip4p_long.cubin $(OBJ_DIR)/lj_tip4p_long.cubin
$(BIN2C) -c -n lj_tip4p_long $(OBJ_DIR)/lj_tip4p_long.cubin > $(OBJ_DIR)/lj_tip4p_long_cubin.h
$(OBJ_DIR)/lal_lj_tip4p_long.o: $(ALL_H) lal_lj_tip4p_long.h lal_lj_tip4p_long.cpp $(OBJ_DIR)/lj_tip4p_long_cubin.h $(OBJ_DIR)/lal_base_atomic.o
$(CUDR) -o $@ -c lal_lj_tip4p_long.cpp -I$(OBJ_DIR)
$(OBJ_DIR)/lal_lj_tip4p_long_ext.o: $(ALL_H) lal_lj_tip4p_long.h lal_lj_tip4p_long_ext.cpp lal_base_atomic.h
$(CUDR) -o $@ -c lal_lj_tip4p_long_ext.cpp -I$(OBJ_DIR)
$(OBJ_DIR)/lj_coul.cubin: lal_lj_coul.cu lal_precision.h lal_preprocessor.h
$(CUDA) --cubin -DNV_KERNEL -o $@ lal_lj_coul.cu

View File

@ -71,7 +71,8 @@ OBJS = $(OBJ_DIR)/lal_atom.o $(OBJ_DIR)/lal_answer.o \
$(OBJ_DIR)/lal_lj_expand_coul_long.o $(OBJ_DIR)/lal_lj_expand_coul_long_ext.o \
$(OBJ_DIR)/lal_coul_long_cs.o $(OBJ_DIR)/lal_coul_long_cs_ext.o \
$(OBJ_DIR)/lal_born_coul_long_cs.o $(OBJ_DIR)/lal_born_coul_long_cs_ext.o \
$(OBJ_DIR)/lal_born_coul_wolf_cs.o $(OBJ_DIR)/lal_born_coul_wolf_cs_ext.o
$(OBJ_DIR)/lal_born_coul_wolf_cs.o $(OBJ_DIR)/lal_born_coul_wolf_cs_ext.o \
$(OBJ_DIR)/lal_lj_tip4p_long.o $(OBJ_DIR)/lal_lj_tip4p_long_ext.o
KERS = $(OBJ_DIR)/device_cl.h $(OBJ_DIR)/atom_cl.h \
$(OBJ_DIR)/neighbor_cpu_cl.h $(OBJ_DIR)/pppm_cl.h \
@ -102,7 +103,8 @@ KERS = $(OBJ_DIR)/device_cl.h $(OBJ_DIR)/atom_cl.h \
$(OBJ_DIR)/lj_cubic_cl.h $(OBJ_DIR)/vashishta_cl.h \
$(OBJ_DIR)/ufm_cl.h $(OBJ_DIR)/dipole_long_lj_cl.h \
$(OBJ_DIR)/lj_expand_coul_long_cl.h $(OBJ_DIR)/coul_long_cs_cl.h \
$(OBJ_DIR)/born_coul_long_cs_cl.h $(OBJ_DIR)/born_coul_wolf_cs_cl.h
$(OBJ_DIR)/born_coul_long_cs_cl.h $(OBJ_DIR)/born_coul_wolf_cs_cl.h \
$(OBJ_DIR)/lj_tip4p_long_cl.h
OCL_EXECS = $(BIN_DIR)/ocl_get_devices
@ -202,6 +204,15 @@ $(OBJ_DIR)/lal_lj.o: $(ALL_H) lal_lj.h lal_lj.cpp $(OBJ_DIR)/lj_cl.h $(OBJ_DIR)
$(OBJ_DIR)/lal_lj_ext.o: $(ALL_H) lal_lj.h lal_lj_ext.cpp lal_base_atomic.h
$(OCL) -o $@ -c lal_lj_ext.cpp -I$(OBJ_DIR)
$(OBJ_DIR)/lj_tip4p_long_cl.h: lal_lj_tip4p_long.cu $(PRE1_H)
$(BSH) ./geryon/file_to_cstr.sh lj_tip4p_long $(PRE1_H) lal_lj_tip4p_long.cu $(OBJ_DIR)/lj_tip4p_long_cl.h;
$(OBJ_DIR)/lal_lj_tip4p_long.o: $(ALL_H) lal_lj_tip4p_long.h lal_lj_tip4p_long.cpp $(OBJ_DIR)/lj_tip4p_long_cl.h $(OBJ_DIR)/lj_tip4p_long_cl.h $(OBJ_DIR)/lal_base_atomic.o
$(OCL) -o $@ -c lal_lj_tip4p_long.cpp -I$(OBJ_DIR)
$(OBJ_DIR)/lal_lj_tip4p_long_ext.o: $(ALL_H) lal_lj_tip4p_long.h lal_lj_tip4p_long_ext.cpp lal_base_atomic.h
$(OCL) -o $@ -c lal_lj_tip4p_long_ext.cpp -I$(OBJ_DIR)
$(OBJ_DIR)/lj_coul_cl.h: lal_lj_coul.cu $(PRE1_H)
$(BSH) ./geryon/file_to_cstr.sh lj_coul $(PRE1_H) lal_lj_coul.cu $(OBJ_DIR)/lj_coul_cl.h;

View File

@ -0,0 +1,241 @@
#if defined(USE_OPENCL)
#include "lj_tip4p_long_cl.h"
#elif defined(USE_CUDART)
const char *lj_tip4p=0;
#else
#include "lj_tip4p_long_cubin.h"
#endif
#include "lal_lj_tip4p_long.h"
#include <cassert>
using namespace LAMMPS_AL;
#define LJ_TIP4PLong_T LJ_TIP4PLong<numtyp, acctyp>
extern Device<PRECISION,ACC_PRECISION> device;
template <class numtyp, class acctyp>
LJ_TIP4PLong<numtyp, acctyp>::LJ_TIP4PLong(): BaseCharge<numtyp,acctyp>(), _allocated(false) {
}
template <class numtyp, class acctyp>
LJ_TIP4PLong<numtyp, acctyp>::~LJ_TIP4PLong() {
clear();
}
template <class numtyp, class acctyp>
int LJ_TIP4PLong<numtyp, acctyp>::bytes_per_atom(const int max_nbors) const {
return this->bytes_per_atom_atomic(max_nbors);
}
template <class numtyp, class acctyp>
int LJ_TIP4PLong<numtyp, acctyp>::init(const int ntypes,
double **host_cutsq, double **host_lj1,
double **host_lj2, double **host_lj3,
double **host_lj4, double **host_offset,
double *host_special_lj, const int nlocal,
const int tH, const int tO,
const double a, const double qd,
const int nall, const int max_nbors,
const int maxspecial, const double cell_size,
const double gpu_split, FILE *_screen,
double **host_cut_ljsq,
const double host_cut_coulsq, const double host_cut_coulsqplus,
double *host_special_coul, const double qqrd2e,
const double g_ewald, int* tag,
int *map_array, int map_size,
int *sametag, int max_same) {
int success;
success=this->init_atomic(nlocal,nall,max_nbors,maxspecial,cell_size,gpu_split,
_screen,lj_tip4p_long,"k_lj_tip4p_long");
if (success!=0)
return success;
k_pair_distrib.set_function(*this->pair_program,"k_lj_tip4p_long_distrib");
TypeH = tH;
TypeO = tO;
alpha = a;
qdist = qd;
// If atom type constants fit in shared memory use fast kernel
int lj_types=ntypes;
shared_types=false;
// int max_shared_types=this->device->max_shared_types();
// if (lj_types<=max_shared_types && this->_block_size>=max_shared_types) {
// lj_types=max_shared_types;
// shared_types=true;
// }
_lj_types=lj_types;
// Allocate a host write buffer for data initialization
UCL_H_Vec<numtyp> host_write(lj_types*lj_types*32,*(this->ucl_device),
UCL_WRITE_ONLY);
for (int i=0; i<lj_types*lj_types; i++)
host_write[i]=0.0;
lj1.alloc(lj_types*lj_types,*(this->ucl_device),UCL_READ_ONLY);
this->atom->type_pack4(ntypes,lj_types,lj1,host_write,host_lj1,host_lj2,
host_cut_ljsq);
lj3.alloc(lj_types*lj_types,*(this->ucl_device),UCL_READ_ONLY);
this->atom->type_pack4(ntypes,lj_types,lj3,host_write,host_lj3,host_lj4,
host_offset);
cutsq.alloc(lj_types*lj_types,*(this->ucl_device),UCL_READ_ONLY);
this->atom->type_pack1(ntypes,lj_types,cutsq,host_write,host_cutsq);
sp_lj.alloc(8,*(this->ucl_device),UCL_READ_ONLY);
for (int i=0; i<4; i++) {
host_write[i]=host_special_lj[i];
host_write[i+4]=host_special_coul[i];
}
ucl_copy(sp_lj,host_write,8,false);
force_comp.alloc(72*72, *(this->ucl_device), UCL_READ_WRITE);
_qqrd2e=qqrd2e;
_g_ewald=g_ewald;
cut_coulsq = host_cut_coulsq;
cut_coulsqplus = host_cut_coulsqplus;
hneight.alloc(nall*4,*(this->ucl_device), UCL_READ_WRITE);
m.alloc(nall,*(this->ucl_device), UCL_READ_WRITE);
ansO.alloc(nall,*(this->ucl_device), UCL_READ_WRITE);
// Allocate a host write buffer for data initialization
UCL_H_Vec<int> host_tag_write(nall,*(this->ucl_device),UCL_WRITE_ONLY);
this->tag.alloc(nall,*(this->ucl_device), UCL_READ_WRITE);
for(int i=0; i<nall; ++i) host_tag_write[i] = tag[i];
ucl_copy(this->tag, host_tag_write, nall, false);
//if(max_same>host_tag_write.cols()) host_tag_write.resize(max_same);
this->atom_sametag.alloc(nall, *(this->ucl_device), UCL_READ_WRITE);
for(int i=0; i<nall; ++i) host_tag_write[i] = sametag[i];
ucl_copy(this->atom_sametag, host_tag_write, nall, false);
if(map_size>host_tag_write.cols()) host_tag_write.resize(map_size);
this->map_array.alloc(map_size,*(this->ucl_device), UCL_READ_WRITE);
for(int i=0; i<map_size; ++i) host_tag_write[i] = map_array[i];
ucl_copy(this->map_array, host_tag_write, map_size, false);
_allocated=true;
this->_max_bytes=lj1.row_bytes()+lj3.row_bytes()+cutsq.row_bytes()+
sp_lj.row_bytes() + hneight.row_bytes()+m.row_bytes()+
this->tag.row_bytes()+this->atom_sametag.row_bytes() +
this->map_array.row_bytes();
return 0;
}
template <class numtyp, class acctyp>
void LJ_TIP4PLong<numtyp, acctyp>::clear() {
if (!_allocated)
return;
_allocated=false;
lj1.clear();
lj3.clear();
sp_lj.clear();
cutsq.clear();
hneight.clear();
m.clear();
tag.clear();
atom_sametag.clear();
map_array.clear();
ansO.clear();
force_comp.clear();
k_pair_distrib.clear();
this->clear_atomic();
}
template <class numtyp, class acctyp>
double LJ_TIP4PLong<numtyp, acctyp>::host_memory_usage() const {
return this->host_memory_usage_atomic()+sizeof(LJ_TIP4PLong<numtyp,acctyp>);
}
// ---------------------------------------------------------------------------
// Calculate energies, forces, and torques
// ---------------------------------------------------------------------------
template <class numtyp, class acctyp>
void LJ_TIP4PLong<numtyp, acctyp>::loop(const bool _eflag, const bool _vflag) {
// Compute the block size and grid size to keep all cores busy
const int BX=this->block_size();
int eflag, vflag;
if (_eflag)
eflag=1;
else
eflag=0;
if (_vflag)
vflag=1;
else
vflag=0;
int GX=static_cast<int>(ceil(static_cast<double>(this->ans->inum())/
(BX/this->_threads_per_atom)));
int ainum=this->ans->inum();
int nbor_pitch=this->nbor->nbor_pitch();
this->time_pair.start();
this->k_pair.set_size(GX,BX);
if (vflag){
this->ansO.resize(ainum*3);
} else {
this->ansO.resize(ainum);
}
this->ansO.zero();
this->k_pair.run(&this->atom->x, &lj1, &lj3, &_lj_types, &sp_lj,
&this->nbor->dev_nbor, &this->_nbor_data->begin(),
&this->ans->force, &this->ans->engv, &eflag, &vflag,
&ainum, &nbor_pitch, &this->_threads_per_atom,
&hneight, &m, &TypeO, &TypeH, &alpha,
&this->atom->q, &cutsq, &_qqrd2e, &_g_ewald,
&cut_coulsq, &cut_coulsqplus, &tag, &map_array,
&atom_sametag, &this->ansO);
GX=static_cast<int>(ceil(static_cast<double>(this->ans->inum())/BX));
this->k_pair_distrib.set_size(GX,BX);
this->k_pair_distrib.run(&this->atom->x, &this->ans->force, &this->ans->engv, &eflag, &vflag,
&ainum, &nbor_pitch, &this->_threads_per_atom,
&hneight, &m, &TypeO, &TypeH, &alpha,
&this->atom->q, &this->ansO);
this->time_pair.stop();
}
template <class numtyp, class acctyp>
void LJ_TIP4PLong<numtyp, acctyp>::copy_relations_data(int **hn, double **newsite, int n,
int* tag, int *map_array, int map_size, int *sametag, int max_same, int ago){
int nall = n;
const int hn_sz = n*4; // matrix size = col size * col number
if(hn_sz > hneight.cols()){
hneight.resize(hn_sz+1);
}
if (ago == 0)
hneight.zero();
if(n > m.cols()){
m.resize(n+1);
}
m.zero();
UCL_H_Vec<int> host_tag_write(nall,*(this->ucl_device),UCL_WRITE_ONLY);
if(this->tag.cols() < nall) this->tag.resize(nall);
for(int i=0; i<nall; ++i) host_tag_write[i] = tag[i];
ucl_copy(this->tag, host_tag_write, nall, false);
if(max_same>host_tag_write.cols()) host_tag_write.resize(max_same);
if(this->atom_sametag.cols() < nall) this->atom_sametag.resize(nall);
for(int i=0; i<nall; ++i) host_tag_write[i] = sametag[i];
ucl_copy(this->atom_sametag, host_tag_write, nall, false);
if(map_size>host_tag_write.cols()) host_tag_write.resize(map_size);
if(this->map_array.cols() < map_size) this->map_array.resize(map_size);
for(int i=0; i<map_size; ++i) host_tag_write[i] = map_array[i];
ucl_copy(this->map_array, host_tag_write, map_size, false);
host_tag_write.clear();
}
template class LJ_TIP4PLong<PRECISION,ACC_PRECISION>;

View File

@ -0,0 +1,518 @@
#ifdef NV_KERNEL
#include "lal_aux_fun1.h"
#ifndef _DOUBLE_DOUBLE
texture<float4> pos_tex;
texture<float> q_tex;
#else
texture<int4,1> pos_tex;
texture<int2> q_tex;
#endif
#else
#define pos_tex x_
#define q_tex q_
#endif
ucl_inline int atom_mapping(const __global int *map, int glob){
return map[glob];
}
ucl_inline int closest_image(int i, int j, const __global int* sametag,
const __global numtyp4 *restrict x_)
{
if (j < 0) return j;
numtyp4 xi; fetch4(xi,i,pos_tex); // = x[i];
numtyp4 xj; fetch4(xj,j,pos_tex);
int closest = j;
numtyp delx = xi.x - xj.x;
numtyp dely = xi.y - xj.y;
numtyp delz = xi.z - xj.z;
numtyp rsqmin = delx*delx + dely*dely + delz*delz;
numtyp rsq;
while (sametag[j] >= 0) {
j = sametag[j];
fetch4(xj,j,pos_tex);
delx = xi.x - xj.x;
dely = xi.y - xj.y;
delz = xi.z - xj.z;
rsq = delx*delx + dely*dely + delz*delz;
if (rsq < rsqmin) {
rsqmin = rsq;
closest = j;
}
}
return closest;
}
ucl_inline void compute_newsite(int iO, int iH1, int iH2,
__global numtyp4 *xM,
numtyp alpha, const __global numtyp4 *restrict x_){
numtyp4 xO; fetch4(xO,iO,pos_tex);
numtyp4 xH1; fetch4(xH1,iH1,pos_tex);
numtyp4 xH2; fetch4(xH2,iH2,pos_tex);
numtyp delx1 = xH1.x - xO.x;
numtyp dely1 = xH1.y - xO.y;
numtyp delz1 = xH1.z - xO.z;
numtyp delx2 = xH2.x - xO.x;
numtyp dely2 = xH2.y - xO.y;
numtyp delz2 = xH2.z - xO.z;
numtyp ap = alpha * (numtyp)0.5;
(*xM).x = xO.x + ap * (delx1 + delx2);
(*xM).y = xO.y + ap * (dely1 + dely2);
(*xM).z = xO.z + ap * (delz1 + delz2);
}
__kernel void k_lj_tip4p_long_distrib(const __global numtyp4 *restrict x_,
__global acctyp4 *restrict ans,
__global acctyp *restrict engv,
const int eflag, const int vflag, const int inum,
const int nbor_pitch, const int t_per_atom,
__global int *restrict hneigh,
__global numtyp4 *restrict m,
const int typeO, const int typeH,
const numtyp alpha,
const __global numtyp *restrict q_, const __global acctyp4 *restrict ansO) {
int tid, ii, offset;
atom_info(t_per_atom,ii,tid,offset);
int i = BLOCK_ID_X*(BLOCK_SIZE_X)+tid;
acctyp4 f;
f.x=(acctyp)0; f.y=(acctyp)0; f.z=(acctyp)0;
if (i<inum) {
numtyp4 ix; fetch4(ix,i,pos_tex);// = x_[i];
int itype = ix.w;
acctyp4 fM, vM;
acctyp eM;
if (itype == typeH) {
int iO = hneigh[i*4];
if (iO < inum) {
fM = ansO[iO];
f.x += fM.x * (acctyp)0.5 * alpha;
f.y += fM.y * (acctyp)0.5 * alpha;
f.z += fM.z * (acctyp)0.5 * alpha;
if (vflag > 0) {
vM = ansO[inum +iO];
engv[inum*2 + i] += vM.x * (acctyp)0.5 * alpha;
engv[inum*3 + i] += vM.y * (acctyp)0.5 * alpha;
engv[inum*4 + i] += vM.z * (acctyp)0.5 * alpha;
vM = ansO[inum*2+iO];
engv[inum*5 + i] += vM.x * (acctyp)0.5 * alpha;
engv[inum*6 + i] += vM.y * (acctyp)0.5 * alpha;
engv[inum*7 + i] += vM.z * (acctyp)0.5 * alpha;
}
}
} else {
fM = ansO[i];
int iH1 = hneigh[i*4 ];
int iH2 = hneigh[i*4+1];
f.x += fM.x * (acctyp)(1 - alpha);
f.y += fM.y * (acctyp)(1 - alpha);
f.z += fM.z * (acctyp)(1 - alpha);
if (eflag > 0) {
eM = engv[i+inum];
engv[inum+i] = eM*(acctyp)(1 - alpha);
if (iH1 < inum) engv[inum+iH1] += eM * (acctyp)0.5 * alpha;
if (iH2 < inum) engv[inum+iH2] += eM * (acctyp)0.5 * alpha;
}
if (vflag > 0) {
vM = ansO[inum + i];
engv[inum*2 + i] += vM.x * (acctyp)(1 - alpha);
engv[inum*3 + i] += vM.y * (acctyp)(1 - alpha);
engv[inum*4 + i] += vM.z * (acctyp)(1 - alpha);
vM = ansO[inum*2 + i];
engv[inum*5 + i] += vM.x * (acctyp)(1 - alpha);
engv[inum*6 + i] += vM.y * (acctyp)(1 - alpha);
engv[inum*7 + i] += vM.z * (acctyp)(1 - alpha);
}
}
acctyp4 old=ans[i];
old.x+=f.x;
old.y+=f.y;
old.z+=f.z;
ans[i]=old;
} // if ii
}
__kernel void k_lj_tip4p_long(const __global numtyp4 *restrict x_,
const __global numtyp4 *restrict lj1,
const __global numtyp4 *restrict lj3,
const int lj_types,
const __global numtyp *restrict sp_lj,
const __global int * dev_nbor,
const __global int * dev_packed,
__global acctyp4 *restrict ans,
__global acctyp *restrict engv,
const int eflag, const int vflag, const int inum,
const int nbor_pitch, const int t_per_atom,
__global int *restrict hneigh,
__global numtyp4 *restrict m,
const int typeO, const int typeH,
const numtyp alpha,
const __global numtyp *restrict q_,
const __global numtyp *restrict cutsq,
const numtyp qqrd2e, const numtyp g_ewald,
const numtyp cut_coulsq, const numtyp cut_coulsqplus,
const __global int *restrict tag, const __global int *restrict map,
const __global int *restrict sametag, __global acctyp4 *restrict ansO) {
int tid, ii, offset;
atom_info(t_per_atom,ii,tid,offset);
acctyp energy = (acctyp)0;
acctyp e_coul = (acctyp)0;
acctyp4 f, fO;
f.x=(acctyp)0; f.y=(acctyp)0; f.z=(acctyp)0;
fO.x=(acctyp)0; fO.y=(acctyp)0; fO.z=(acctyp)0;
acctyp virial[6],vO[6];
for (int i=0; i<6; i++) {
virial[i]=(acctyp)0;
vO[i]=(acctyp)0;
}
if (ii<inum) {
int i, numj, nbor, nbor_end;
__local int n_stride;
nbor_info(dev_nbor,dev_packed,nbor_pitch,t_per_atom,ii,offset,i,numj,
n_stride,nbor_end,nbor);
numtyp4 ix; fetch4(ix,i,pos_tex); //x_[i];
numtyp qtmp; fetch(qtmp,i,q_tex);
int itype = ix.w;
numtyp4 x1 = ix;
int non_local_oxy = 0;
int iH1, iH2, iO;
if(itype == typeO) {
iO = i;
if (hneigh[i*4+2] != -1) {
iH1 = atom_mapping(map, tag[i] + 1);
iH2 = atom_mapping(map, tag[i] + 2);
// set iH1,iH2 to closest image to O
iH1 = closest_image(i, iH1, sametag, x_);
iH2 = closest_image(i, iH2, sametag, x_);
hneigh[ i*4 ] = iH1;
hneigh[ i*4+1] = iH2;
hneigh[ i*4+2] = -1;
hneigh[iH1*4 ] = i;
hneigh[iH1*4+1] += -1;
hneigh[iH1*4+2] = -1;
hneigh[iH2*4 ] = i;
hneigh[iH2*4+1] += -1;
hneigh[iH2*4+2] = -1;
} else {
iH1 = hneigh[i*4 ];
iH2 = hneigh[i*4+1];
}
if(m[iO].w == 0) {
compute_newsite(iO,iH1,iH2, &m[iO], alpha, x_);
m[iO].w = qtmp;
}
x1 = m[iO];
} else {
if (hneigh[i*4+2] != -1) {
int iI, iH;
iI = atom_mapping(map,tag[i] - 1);
numtyp4 iIx; fetch4(iIx,iI,pos_tex); //x_[iI];
if (iIx.w == typeH) {
iO = atom_mapping(map,tag[i] - 2);
iO = closest_image(i, iO, sametag, x_);
iH1 = closest_image(i, iI, sametag, x_);
iH2 = i;
} else if (iIx.w == typeO) {
iH = atom_mapping(map, tag[i] + 1);
iO = closest_image(i,iI,sametag, x_);
iH1 = i;
iH2 = closest_image(i,iH,sametag, x_);
}
hneigh[iH1*4+0] = iO;
hneigh[iH1*4+1] += -1;
hneigh[iH1*4+2] = -1;
hneigh[iH2*4+0] = iO;
hneigh[iH2*4+1] += -1;
hneigh[iH2*4+2] = -1;
hneigh[iO*4+0] = iH1;
hneigh[iO*4+1] = iH2;
hneigh[iO*4+2] = -1;
} else {
iO = hneigh[i *4 ];
iH1 = hneigh[iO*4 ];
iH2 = hneigh[iO*4+1];
}
if (iO >= inum) {
non_local_oxy = 1;
if(m[iO].w == 0) {
compute_newsite(iO,iH1,iH2, &m[iO], alpha, x_);
numtyp qO; fetch(qO,iO,q_tex);
m[iO].w = qO;
}
}
}
for ( ; nbor<nbor_end; nbor+=n_stride) {
int j=dev_packed[nbor];
numtyp factor_lj,factor_coul;
factor_lj = sp_lj[sbmask(j)];
factor_coul = (numtyp)1.0-sp_lj[sbmask(j)+4];
j &= NEIGHMASK;
numtyp4 jx; fetch4(jx,j,pos_tex); //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 rsq = delx*delx+dely*dely+delz*delz;
int mtype = itype*lj_types+jtype;
if (rsq < lj1[mtype].z) { // cut_ljsq
numtyp r2inv = ucl_recip(rsq);
numtyp r6inv = r2inv*r2inv*r2inv;
numtyp forcelj = r6inv*(lj1[mtype].x*r6inv-lj1[mtype].y);
forcelj *= r2inv*factor_lj;
f.x += delx*forcelj;
f.y += dely*forcelj;
f.z += delz*forcelj;
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*forcelj;
virial[1] += dely*dely*forcelj;
virial[2] += delz*delz*forcelj;
virial[3] += delx*dely*forcelj;
virial[4] += delx*delz*forcelj;
virial[5] += dely*delz*forcelj;
}
} // if LJ
if (rsq < cut_coulsqplus) { //cut_coulsqplus
int jH1, jH2, jO;
numtyp qj; fetch(qj,j,q_tex);
numtyp4 x2 = jx;
if(itype == typeO || jtype == typeO) {
if (jtype == typeO) {
jO = j;
if (hneigh[j*4+2] != -1) {
jH1 = atom_mapping(map,tag[j] + 1);
jH2 = atom_mapping(map,tag[j] + 2);
// set iH1,iH2 to closest image to O
jH1 = closest_image(j, jH1, sametag, x_);
jH2 = closest_image(j, jH2, sametag, x_);
hneigh[j*4 ] = jH1;
hneigh[j*4+1] = jH2;
hneigh[j*4+2] = -1;
hneigh[jH1*4 ] = j;
hneigh[jH1*4+1] += -1;
hneigh[jH1*4+2] = -1;
hneigh[jH2*4 ] = j;
hneigh[jH2*4+1] += -1;
hneigh[jH2*4+2] = -1;
} else {
jH1 = hneigh[j*4 ];
jH2 = hneigh[j*4+1];
}
if (m[j].w == 0) {
compute_newsite(j, jH1, jH2, &m[j], alpha, x_);
m[j].w = qj;
}
x2 = m[j];
}
delx = x1.x-x2.x;
dely = x1.y-x2.y;
delz = x1.z-x2.z;
rsq = delx*delx+dely*dely+delz*delz;
}
if (rsq < cut_coulsq) {
numtyp r2inv = ucl_recip(rsq);
numtyp r = ucl_rsqrt(r2inv);
numtyp grij = g_ewald * r;
numtyp expm2 = ucl_exp(-grij*grij);
numtyp t = ucl_recip((numtyp)1.0 + EWALD_P*grij);
numtyp _erfc = t * (A1+t*(A2+t*(A3+t*(A4+t*A5)))) * expm2;
numtyp prefactor = qj;
prefactor *= qqrd2e*qtmp/r;
numtyp force_coul = r2inv*prefactor * (_erfc + EWALD_F*grij*expm2 - factor_coul);
if (itype == typeH) {
f.x += delx * force_coul;
f.y += dely * force_coul;
f.z += delz * force_coul;
} else {
fO.x += delx * force_coul;
fO.y += dely * force_coul;
fO.z += delz * force_coul;
fO.w += -1;
}
if (eflag>0) {
e_coul += prefactor*(_erfc-factor_coul);
}
if (vflag>0) {
acctyp4 fd;
fd.x = delx*force_coul;
fd.y = dely*force_coul;
fd.z = delz*force_coul;
if (itype == typeH) {
if (jtype == typeH){
virial[0] += delx*fd.x;
virial[1] += dely*fd.y;
virial[2] += delz*fd.z;
virial[3] += delx*fd.y;
virial[4] += delx*fd.z;
virial[5] += dely*fd.z;
} else {
numtyp cO = 1 - alpha, cH = 0.5*alpha;
numtyp4 vdj;
numtyp4 xjH1; fetch4(xjH1,jH1,pos_tex);
numtyp4 xjH2; fetch4(xjH2,jH2,pos_tex);
numtyp4 xjO; fetch4(xjO,jO,pos_tex);
vdj.x = xjO.x*cO + xjH1.x*cH + xjH2.x*cH;
vdj.y = xjO.y*cO + xjH1.y*cH + xjH2.y*cH;
vdj.z = xjO.z*cO + xjH1.z*cH + xjH2.z*cH;
virial[0] += (ix.x - vdj.x)*fd.x;
virial[1] += (ix.y - vdj.y)*fd.y;
virial[2] += (ix.z - vdj.z)*fd.z;
virial[3] += (ix.x - vdj.x)*fd.y;
virial[4] += (ix.x - vdj.x)*fd.z;
virial[5] += (ix.y - vdj.y)*fd.z;
}
} else {
numtyp cO = 1 - alpha, cH = 0.5*alpha;
numtyp4 vdi, vdj;
numtyp4 xH1; fetch4(xH1,iH1,pos_tex);
numtyp4 xH2; fetch4(xH2,iH2,pos_tex);
numtyp4 xO; fetch4(xO,iO,pos_tex);
vdi.x = xO.x*cO + xH1.x*cH + xH2.x*cH;
vdi.y = xO.y*cO + xH1.y*cH + xH2.y*cH;
vdi.z = xO.z*cO + xH1.z*cH + xH2.z*cH;
if (jtype != typeH){
numtyp4 xjH1; fetch4(xjH1,jH1,pos_tex);
numtyp4 xjH2; fetch4(xjH2,jH2,pos_tex);
numtyp4 xjO; fetch4(xjO,jO,pos_tex);
vdj.x = xjO.x*cO + xjH1.x*cH + xjH2.x*cH;
vdj.y = xjO.y*cO + xjH1.y*cH + xjH2.y*cH;
vdj.z = xjO.z*cO + xjH1.z*cH + xjH2.z*cH;
} else vdj = jx;
vO[0] += 0.5*(vdi.x - vdj.x)*fd.x;
vO[1] += 0.5*(vdi.y - vdj.y)*fd.y;
vO[2] += 0.5*(vdi.z - vdj.z)*fd.z;
vO[3] += 0.5*(vdi.x - vdj.x)*fd.y;
vO[4] += 0.5*(vdi.x - vdj.x)*fd.z;
vO[5] += 0.5*(vdi.y - vdj.y)*fd.z;
}
}
}
if (non_local_oxy == 1) {
if (iO == j) {
x2 = ix;
qj = qtmp;
}
numtyp4 x1m = m[iO];
delx = x1m.x-x2.x;
dely = x1m.y-x2.y;
delz = x1m.z-x2.z;
rsq = delx*delx+dely*dely+delz*delz;
if (rsq < cut_coulsq) {
numtyp r2inv = ucl_recip(rsq);
numtyp r = ucl_rsqrt(r2inv);
numtyp grij = g_ewald * r;
numtyp expm2 = ucl_exp(-grij*grij);
numtyp t = ucl_recip((numtyp)1.0 + EWALD_P*grij);
numtyp _erfc = t * (A1+t*(A2+t*(A3+t*(A4+t*A5)))) * expm2;
numtyp prefactor = qj;
prefactor *= qqrd2e*x1m.w/r;
numtyp force_coul = r2inv*prefactor * (_erfc + EWALD_F*grij*expm2 - factor_coul);
numtyp cO = 1 - alpha, cH = 0.5*alpha;
numtyp4 fd;
fd.x = delx * force_coul * cH;
fd.y = dely * force_coul * cH;
fd.z = delz * force_coul * cH;
f.x += fd.x;
f.y += fd.y;
f.z += fd.z;
if (eflag>0) {
e_coul += prefactor*(_erfc-factor_coul) * (acctyp)0.5 * alpha;
}
if (vflag>0) {
numtyp4 xH1; fetch4(xH1,iH1,pos_tex);
numtyp4 xH2; fetch4(xH2,iH2,pos_tex);
numtyp4 xO; fetch4(xO,iO,pos_tex);
virial[0] += ((xO.x*cO + xH1.x*cH + xH2.x*cH) - x2.x) * fd.x;
virial[1] += ((xO.y*cO + xH1.y*cH + xH2.y*cH) - x2.y) * fd.y;
virial[2] += ((xO.z*cO + xH1.z*cH + xH2.z*cH) - x2.z) * fd.z;
virial[3] += ((xO.x*cO + xH1.x*cH + xH2.x*cH) - x2.x) * fd.y;
virial[4] += ((xO.x*cO + xH1.x*cH + xH2.x*cH) - x2.x) * fd.z;
virial[5] += ((xO.y*cO + xH1.y*cH + xH2.y*cH) - x2.y) * fd.z;
}
}
}
} // if cut_coulsqplus
} // for nbor
if (t_per_atom>1) {
__local acctyp red_acc[6][BLOCK_PAIR];
red_acc[0][tid]=fO.x;
red_acc[1][tid]=fO.y;
red_acc[2][tid]=fO.z;
red_acc[3][tid]=fO.w;
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];
}
}
fO.x=red_acc[0][tid];
fO.y=red_acc[1][tid];
fO.z=red_acc[2][tid];
fO.w=red_acc[3][tid];
if (vflag>0) {
for (int r=0; r<6; r++) red_acc[r][tid]=vO[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++) vO[r]=red_acc[r][tid];
}
}
if(offset == 0) {
ansO[i] = fO;
if (vflag>0) {
ansO[inum + i].x = vO[0];
ansO[inum + i].y = vO[1];
ansO[inum + i].z = vO[2];
ansO[inum*2 + i].x = vO[3];
ansO[inum*2 + i].y = vO[4];
ansO[inum*2 + i].z = vO[5];
}
}
store_answers_q(f,energy,e_coul,virial,ii,inum,tid,t_per_atom,offset,eflag,
vflag,ans,engv);
} // if ii
}
__kernel void k_lj_tip4p_long_fast(){}

View File

@ -0,0 +1,95 @@
#ifndef LAL_LJ_TIP4P_LONG_H
#define LAL_LJ_TIP4P_LONG_H
#include "lal_base_charge.h"
namespace LAMMPS_AL {
template <class numtyp, class acctyp>
class LJ_TIP4PLong : public BaseCharge<numtyp, acctyp> {
public:
LJ_TIP4PLong();
~LJ_TIP4PLong();
/// Clear any previous data and set up for a new LAMMPS run
/** \param max_nbors initial number of rows in the neighbor matrix
* \param cell_size cutoff + skin
* \param gpu_split fraction of particles handled by device
*
* Returns:
* - 0 if successfull
* - -1 if fix gpu not found
* - -3 if there is an out of memory error
* - -4 if the GPU library was not compiled for GPU
* - -5 Double precision is not supported on card **/
int init(const int ntypes, double **host_cutsq,
double **host_lj1, double **host_lj2, double **host_lj3,
double **host_lj4, double **host_offset, double *host_special_lj,
const int nlocal, const int tH, const int tO,
const double alpha, const double qdist,
const int nall, const int max_nbors,
const int maxspecial, const double cell_size,
const double gpu_split, FILE *screen,
double **host_cut_ljsq,
const double host_cut_coulsq, const double host_cut_coulsqplus,
double *host_special_coul, const double qqrd2e,
const double g_ewald, int* tag,
int *map_array, int map_size,
int *sametag, int max_same);
/// Clear all host and device data
/** \note This is called at the beginning of the init() routine **/
void clear();
/// Returns memory usage on device per atom
int bytes_per_atom(const int max_nbors) const;
/// Total host memory used by library for pair style
double host_memory_usage() const;
/// Copy data which is easier to compute in LAMMPS_NS
void copy_relations_data(int **hn, double **m, int n,int* tag, int *map_array, int map_size,
int *sametag, int max_same, int ago);
// --------------------------- TYPE DATA --------------------------
/// lj1.x = lj1, lj1.y = lj2, lj1.z = cutsq_vdw
UCL_D_Vec<numtyp4> lj1;
/// lj3.x = lj3, lj3.y = lj4, lj3.z = offset
UCL_D_Vec<numtyp4> lj3;
/// cutsq
UCL_D_Vec<numtyp> cutsq;
/// Special LJ values [0-3] and Special Coul values [4-7]
UCL_D_Vec<numtyp> sp_lj;
/// If atom type constants fit in shared memory, use fast kernels
bool shared_types;
/// Number of atom types
int _lj_types;
numtyp _qqrd2e, _g_ewald;
/// TIP4P water parameters
int TypeO, TypeH;
numtyp alpha, qdist;
numtyp cut_coulsq, cut_coulsqplus;
UCL_D_Vec<int> hneight;
UCL_D_Vec<numtyp4> m; // position and charge of virtual particle
UCL_D_Vec<acctyp4> ansO; // force applied to virtual particle
UCL_D_Vec<acctyp4> force_comp;
UCL_D_Vec<int> tag;
UCL_D_Vec<int> map_array;
UCL_D_Vec<int> atom_sametag;
UCL_Kernel k_pair_distrib;
private:
bool _allocated;
void loop(const bool _eflag, const bool _vflag);
};
}
#endif

View File

@ -0,0 +1,135 @@
#include <iostream>
#include <cassert>
#include <cmath>
#include "lal_lj_tip4p_long.h"
using namespace std;
using namespace LAMMPS_AL;
static LJ_TIP4PLong<PRECISION,ACC_PRECISION> LJTIP4PLMF;
// ---------------------------------------------------------------------------
// Allocate memory on host and device and copy constants to device
// ---------------------------------------------------------------------------
int ljtip4p_long_gpu_init(const int ntypes, double **cutsq, double **host_lj1,
double **host_lj2, double **host_lj3, double **host_lj4,
double **offset, double *special_lj, const int inum,
const int tH, const int tO,
const double alpha, const double qdist,
const int nall, const int max_nbors, const int maxspecial,
const double cell_size, int &gpu_mode, FILE *screen,
double **host_cut_ljsq,
const double host_cut_coulsq, const double host_cut_coulsqplus,
double *host_special_coul, const double qqrd2e,
const double g_ewald, int* tag,
int *map_array, int map_size,
int *sametag, int max_same) {
LJTIP4PLMF.clear();
gpu_mode=LJTIP4PLMF.device->gpu_mode();
double gpu_split=LJTIP4PLMF.device->particle_split();
int first_gpu=LJTIP4PLMF.device->first_device();
int last_gpu=LJTIP4PLMF.device->last_device();
int world_me=LJTIP4PLMF.device->world_me();
int gpu_rank=LJTIP4PLMF.device->gpu_rank();
int procs_per_gpu=LJTIP4PLMF.device->procs_per_gpu();
LJTIP4PLMF.device->init_message(screen,"lj/cut/tip4p/long/gpu",first_gpu,last_gpu);
bool message=false;
if (LJTIP4PLMF.device->replica_me()==0 && screen)
message=true;
if (message) {
fprintf(screen,"Initializing Device and compiling on process 0...");
fflush(screen);
}
int init_ok=0;
if (world_me==0)
init_ok=LJTIP4PLMF.init(ntypes, cutsq, host_lj1, host_lj2, host_lj3,
host_lj4, offset, special_lj, inum,
tH, tO, alpha, qdist, nall, 300,
maxspecial, cell_size, gpu_split, screen,
host_cut_ljsq, host_cut_coulsq, host_cut_coulsqplus,
host_special_coul, qqrd2e, g_ewald, tag,
map_array, map_size,
sametag, max_same);
LJTIP4PLMF.device->world_barrier();
if (message)
fprintf(screen,"Done.\n");
for (int i=0; i<procs_per_gpu; i++) {
if (message) {
if (last_gpu-first_gpu==0)
fprintf(screen,"Initializing Device %d on core %d...",first_gpu,i);
else
fprintf(screen,"Initializing Devices %d-%d on core %d...",first_gpu,
last_gpu,i);
fflush(screen);
}
if (gpu_rank==i && world_me!=0)
init_ok=LJTIP4PLMF.init(ntypes, cutsq, host_lj1, host_lj2, host_lj3, host_lj4,
offset, special_lj, inum,
tH, tO, alpha, qdist, nall, 300, maxspecial,
cell_size, gpu_split, screen, host_cut_ljsq,
host_cut_coulsq, host_cut_coulsqplus,
host_special_coul, qqrd2e, g_ewald,tag,
map_array, map_size,
sametag, max_same);
LJTIP4PLMF.device->gpu_barrier();
if (message)
fprintf(screen,"Done.\n");
}
if (message)
fprintf(screen,"\n");
if (init_ok==0)
LJTIP4PLMF.estimate_gpu_overhead();
return init_ok;
}
void ljtip4p_long_gpu_clear() {
LJTIP4PLMF.clear();
}
int ** ljtip4p_long_gpu_compute_n(const int ago, const int inum_full,
const int nall, double **host_x, int *host_type,
double *sublo, double *subhi, tagint *tag, int **nspecial,
tagint **special, const bool eflag, const bool vflag,
const bool eatom, const bool vatom, int &host_start,
int **ilist, int **jnum, const double cpu_time,
bool &success, double *host_q, double *boxlo,
double *prd) {
return LJTIP4PLMF.compute(ago, inum_full, nall, host_x, host_type, sublo,
subhi, tag, nspecial, special, eflag, vflag, eatom,
vatom, host_start, ilist, jnum, cpu_time, success,
host_q,boxlo, prd);
}
void ljtip4p_long_gpu_compute(const int ago, const int inum_full, const int nall,
double **host_x, int *host_type, int *ilist, int *numj,
int **firstneigh, const bool eflag, const bool vflag,
const bool eatom, const bool vatom, int &host_start,
const double cpu_time, bool &success,double *host_q,
const int nlocal, double *boxlo, double *prd) {
LJTIP4PLMF.compute(ago,inum_full,nall,host_x,host_type,ilist,numj,
firstneigh,eflag,vflag,eatom,vatom,host_start,cpu_time,success,host_q,
nlocal,boxlo,prd);
}
double ljtip4p_long_gpu_bytes() {
return LJTIP4PLMF.host_memory_usage();
}
void ljtip4p_long_copy_molecule_data(int **hn, double **m, int n, int* tag,
int *map_array, int map_size,
int *sametag, int max_same, int ago){
LJTIP4PLMF.copy_relations_data(hn, m, n, tag,map_array,map_size,sametag, max_same, ago);
}

View File

@ -143,6 +143,8 @@ action pair_ufm_gpu.cpp
action pair_ufm_gpu.h
action pair_lj_cut_dipole_long_gpu.cpp pair_lj_cut_dipole_long.cpp
action pair_lj_cut_dipole_long_gpu.h pair_lj_cut_dipole_long.cpp
action pair_lj_cut_tip4p_long_gpu.h
action pair_lj_cut_tip4p_long_gpu.cpp
# edit 2 Makefile.package files to include/exclude package info

View File

@ -0,0 +1,234 @@
#include <cmath>
#include <cstdio>
#include <cstdlib>
#include <cstring>
#include "pair_lj_cut_tip4p_long_gpu.h"
#include "atom.h"
#include "atom_vec.h"
#include "comm.h"
#include "force.h"
#include "neighbor.h"
#include "neigh_list.h"
#include "integrate.h"
#include "memory.h"
#include "error.h"
#include "neigh_request.h"
#include "universe.h"
#include "update.h"
#include "domain.h"
#include "kspace.h"
#include "angle.h"
#include "bond.h"
#include "gpu_extra.h"
#define EWALD_F 1.12837917
#define EWALD_P 0.3275911
#define A1 0.254829592
#define A2 -0.284496736
#define A3 1.421413741
#define A4 -1.453152027
#define A5 1.061405429
using namespace LAMMPS_NS;
// External functions from cuda library for atom decomposition
int ljtip4p_long_gpu_init(const int ntypes, double **cutsq, double **host_lj1,
double **host_lj2, double **host_lj3, double **host_lj4,
double **offset, double *special_lj, const int nlocal,
const int tH, const int tO, const double alpha, const double qdist,
const int nall, const int max_nbors, const int maxspecial,
const double cell_size, int &gpu_mode, FILE *screen,
double **host_cut_ljsq, const double host_cut_coulsq,
const double host_cut_coulsqplus,
double *host_special_coul, const double qqrd2e,
const double g_ewald, int* tag,
int *map_array, int map_size,
int *sametag, int max_same);
void ljtip4p_long_gpu_clear();
int ** ljtip4p_long_gpu_compute_n(const int ago, const int inum,
const int nall, double **host_x, int *host_type,
double *sublo, double *subhi, tagint *tag, int **nspecial,
tagint **special, const bool eflag, const bool vflag,
const bool eatom, const bool vatom, int &host_start,
int **ilist, int **jnum,
const double cpu_time, bool &success, double *host_q,
double *boxlo, double *prd);
void ljtip4p_long_gpu_compute(const int ago, const int inum, const int nall,
double **host_x, int *host_type, int *ilist, int *numj,
int **firstneigh, const bool eflag, const bool vflag,
const bool eatom, const bool vatom, int &host_start,
const double cpu_time,
bool &success, double *host_q, const int nlocal,
double *boxlo, double *prd);
double ljtip4p_long_gpu_bytes();
void ljtip4p_long_copy_molecule_data(int **, double **, int, int* , int *, int, int *, int , int);
/* ---------------------------------------------------------------------- */
PairLJCutTIP4PLongGPU::PairLJCutTIP4PLongGPU(LAMMPS *lmp)
: PairLJCutTIP4PLong(lmp), gpu_mode(GPU_FORCE)
{
respa_enable = 0;
reinitflag = 0;
cpu_time = 0.0;
GPU_EXTRA::gpu_ready(lmp->modify, lmp->error);
}
/* ----------------------------------------------------------------------
free all arrays
------------------------------------------------------------------------- */
PairLJCutTIP4PLongGPU::~PairLJCutTIP4PLongGPU()
{
ljtip4p_long_gpu_clear();
}
/* ---------------------------------------------------------------------- */
void PairLJCutTIP4PLongGPU::compute(int eflag, int vflag)
{
if (eflag || vflag) ev_setup(eflag,vflag);
else evflag = vflag_fdotr = 0;
int nall = atom->nlocal + atom->nghost;
int inum, host_start;
ljtip4p_long_copy_molecule_data(hneigh, newsite, nall, atom->tag,
atom->get_map_array(), atom->get_map_size(),
atom->sametag, atom->get_max_same(), neighbor->ago);
bool success = true;
int *ilist, *numneigh, **firstneigh;
if (gpu_mode != GPU_FORCE) {
inum = atom->nlocal;
firstneigh = ljtip4p_long_gpu_compute_n(neighbor->ago, inum, nall,
atom->x, atom->type, domain->sublo,
domain->subhi, atom->tag, atom->nspecial,
atom->special, eflag, vflag, eflag_atom,
vflag_atom, host_start, &ilist, &numneigh,
cpu_time, success, atom->q, domain->boxlo,
domain->prd);
} else {
inum = list->inum;
ilist = list->ilist;
numneigh = list->numneigh;
firstneigh = list->firstneigh;
ljtip4p_long_gpu_compute(neighbor->ago, inum, nall, atom->x, atom->type,
ilist, numneigh, firstneigh, eflag, vflag, eflag_atom,
vflag_atom, host_start, cpu_time, success, atom->q,
atom->nlocal, domain->boxlo, domain->prd);
}
if (!success)
error->one(FLERR,"Insufficient memory on accelerator");
// if (host_start<inum) {
// cpu_time = MPI_Wtime();
// cpu_compute(host_start, inum, eflag, vflag, ilist, numneigh, firstneigh);
// cpu_time = MPI_Wtime() - cpu_time;
// }
}
/* ----------------------------------------------------------------------
init specific to this pair style
------------------------------------------------------------------------- */
void PairLJCutTIP4PLongGPU::init_style()
{
cut_respa = NULL;
if (atom->tag_enable == 0)
error->all(FLERR,"Pair style lj/cut/tip4p/long/gpu requires atom IDs");
if (!atom->q_flag)
error->all(FLERR,
"Pair style lj/cut/tip4p/long/gpu requires atom attribute q");
if (force->bond == NULL)
error->all(FLERR,"Must use a bond style with TIP4P potential");
if (force->angle == NULL)
error->all(FLERR,"Must use an angle style with TIP4P potential");
if (atom->map_style == 2)
error->all(FLERR,"GPU-accelerated lj/cut/tip4p/long currently requires map style 'array' (atom_modify map array)");
//PairLJCutCoulLong::init_style();
// Repeat cutsq calculation because done after call to init_style
double maxcut = -1.0;
double cut;
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);
cut *= cut;
if (cut > maxcut)
maxcut = cut;
cutsq[i][j] = cutsq[j][i] = cut;
} else
cutsq[i][j] = cutsq[j][i] = 0.0;
}
}
double cell_size = sqrt(maxcut) + neighbor->skin;
// insure use of KSpace long-range solver, set g_ewald
if (force->kspace == NULL)
error->all(FLERR,"Pair style requires a KSpace style");
g_ewald = force->kspace->g_ewald;
// setup force tables
if (ncoultablebits) init_tables(cut_coul,cut_respa);
int maxspecial=0;
if (atom->molecular)
maxspecial=atom->maxspecial;
// set alpha parameter
double theta = force->angle->equilibrium_angle(typeA);
double blen = force->bond->equilibrium_distance(typeB);
alpha = qdist / (cos(0.5*theta) * blen);
cut_coulsq = cut_coul * cut_coul;
double cut_coulsqplus = (cut_coul+qdist+blen) * (cut_coul+qdist+blen);
if (maxcut < cut_coulsqplus) {
cell_size = (cut_coul+qdist+blen) + neighbor->skin;
}
if (comm->cutghostuser < cell_size) {
comm->cutghostuser = cell_size;
if (comm->me == 0)
error->warning(FLERR,"Increasing communication cutoff for TIP4P GPU style");
}
int success = ljtip4p_long_gpu_init(atom->ntypes+1, cutsq, lj1, lj2, lj3, lj4,
offset, force->special_lj, atom->nlocal,
typeH, typeO, alpha, qdist,
atom->nlocal+atom->nghost, 300, maxspecial,
cell_size, gpu_mode, screen, cut_ljsq,
cut_coulsq, cut_coulsqplus,
force->special_coul, force->qqrd2e,
g_ewald,
atom->tag,
atom->get_map_array(), atom->get_map_size(),
atom->sametag, atom->get_max_same());
GPU_EXTRA::check_flag(success,error,world);
if (gpu_mode == GPU_FORCE) {
int irequest = neighbor->request(this,instance_me);
neighbor->requests[irequest]->half = 0;
neighbor->requests[irequest]->full = 1;
neighbor->requests[irequest]->cut = 1;
neighbor->requests[irequest]->cutoff = cut_coul+qdist+blen + neighbor->skin;
}
}
/* ---------------------------------------------------------------------- */
double PairLJCutTIP4PLongGPU::memory_usage()
{
double bytes = PairLJCutTIP4PLong::memory_usage();
return bytes + ljtip4p_long_gpu_bytes();
}
/* ---------------------------------------------------------------------- */
//void PairLJCutTIP4PLongGPU::cpu_compute(int start, int inum, int eflag, int vflag,
// int *ilist, int *numneigh, int **firstneigh) {
// error->all(FLERR,"PairLJCutTIP4PLongGPU::cpu_compute not implemented");
//}

View File

@ -0,0 +1,32 @@
#ifdef PAIR_CLASS
PairStyle(lj/cut/tip4p/long/gpu,PairLJCutTIP4PLongGPU)
#else
#ifndef LMP_PAIR_LJ_TIP4P_LONG_GPU_H
#define LMP_PAIR_LJ_TIP4P_LONG_GPU_H
#include "pair_lj_cut_tip4p_long.h"
namespace LAMMPS_NS {
class PairLJCutTIP4PLongGPU : public PairLJCutTIP4PLong {
public:
PairLJCutTIP4PLongGPU(LAMMPS *lmp);
~PairLJCutTIP4PLongGPU();
// void cpu_compute(int, int, int, int, int *, int *, int **);
void compute(int, int);
void init_style();
double memory_usage();
enum { GPU_FORCE, GPU_NEIGH, GPU_HYB_NEIGH };
private:
int gpu_mode;
double cpu_time;
};
}
#endif
#endif

View File

@ -287,6 +287,7 @@ class Atom : protected Pointers {
inline int* get_map_array() {return map_array;};
inline int get_map_size() {return map_tag_max+1;};
inline int get_max_same() {return max_same;};
inline int get_map_maxarray() {return map_maxarray+1;};
bigint memory_usage();