1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402
|
/**************************************************************************
lj_tip4p_long.cpp
-------------------
V. Nikolskiy (HSE)
Class for acceleration of the lj/tip4p/long pair style
__________________________________________________________________________
This file is part of the LAMMPS Accelerator Library (LAMMPS_AL)
__________________________________________________________________________
begin :
email : thevsevak@gmail.com
***************************************************************************/
#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>
namespace LAMMPS_AL {
#define LJTIP4PLongT LJ_TIP4PLong<numtyp, acctyp>
extern Device<PRECISION,ACC_PRECISION> device;
template <class numtyp, class acctyp>
LJTIP4PLongT::LJ_TIP4PLong(): BaseCharge<numtyp,acctyp>(), _allocated(false) {
}
template <class numtyp, class acctyp>
LJTIP4PLongT::~LJ_TIP4PLong() {
clear();
}
template <class numtyp, class acctyp>
int LJTIP4PLongT::bytes_per_atom(const int max_nbors) const {
return this->bytes_per_atom_atomic(max_nbors);
}
template <class numtyp, class acctyp>
int LJTIP4PLongT::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 map_size, 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");
k_pair_reneigh.set_function(*this->pair_program,"k_lj_tip4p_reneigh");
k_pair_newsite.set_function(*this->pair_program,"k_lj_tip4p_newsite");
#if defined(LAL_OCL_EV_JIT)
k_pair_distrib_noev.set_function(*this->pair_program_noev,
"k_lj_tip4p_long_distrib");
#else
k_pair_dt_sel = &k_pair_distrib;
#endif
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);
this->tag.alloc(nall,*(this->ucl_device), UCL_READ_ONLY);
this->atom_sametag.alloc(max_same, *(this->ucl_device), UCL_READ_ONLY);
this->map_array.alloc(map_size,*(this->ucl_device), UCL_READ_ONLY);
_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 LJTIP4PLongT::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();
k_pair_reneigh.clear();
k_pair_newsite.clear();
#if defined(LAL_OCL_EV_JIT)
k_pair_distrib_noev.clear();
#endif
this->clear_atomic();
}
template <class numtyp, class acctyp>
double LJTIP4PLongT::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>
int LJTIP4PLongT::loop(const int eflag, const int vflag) {
// Compute the block size and grid size to keep all cores busy
const int BX=this->block_size();
int ainum=this->ans->inum();
const int nall = this->atom->nall();
int nbor_pitch=this->nbor->nbor_pitch();
this->time_pair.start();
int GX;
GX=static_cast<int>(ceil(static_cast<double>(nall)/BX));
if (t_ago == 0) {
this->k_pair_reneigh.set_size(GX,BX);
this->k_pair_reneigh.run(&this->atom->x,
&this->nbor->dev_nbor, &this->_nbor_data->begin(),
&nall, &ainum,&nbor_pitch, &this->_threads_per_atom,
&hneight, &m, &TypeO, &TypeH,
&tag, &map_array, &atom_sametag);
}
this->k_pair_newsite.set_size(GX,BX);
this->k_pair_newsite.run(&this->atom->x,
&this->nbor->dev_nbor, &this->_nbor_data->begin(),
&nall, &ainum,
&nbor_pitch, &this->_threads_per_atom,
&hneight, &m, &TypeO, &TypeH, &alpha,
&this->atom->q);
GX=static_cast<int>(ceil(static_cast<double>(this->ans->inum())/
(BX/this->_threads_per_atom)));
if (vflag) {
this->ansO.resize_ib(ainum*3);
} else {
this->ansO.resize_ib(ainum);
}
this->ansO.zero();
this->device->gpu->sync();
if(shared_types) {
this->k_pair_sel->set_size(GX,BX);
this->k_pair_sel->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, &this->ansO);
} else {
this->k_pair.set_size(GX,BX);
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, &this->ansO);
}
#if defined(LAL_OCL_EV_JIT)
if (eflag || vflag) k_pair_dt_sel = &k_pair_distrib;
else k_pair_dt_sel = &k_pair_distrib_noev;
#endif
GX=static_cast<int>(ceil(static_cast<double>(this->ans->inum())/BX));
k_pair_dt_sel->set_size(GX,BX);
k_pair_dt_sel->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();
return GX;
}
template <class numtyp, class acctyp>
void LJTIP4PLongT::copy_relations_data(int n, tagint *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
hneight.resize_ib(hn_sz);
m.resize_ib(n);
m.zero();
if (ago == 0) {
hneight.zero();
{
UCL_H_Vec<tagint> host_tag_write;
host_tag_write.view(tag, nall, *(this->ucl_device));
this->tag.resize_ib(nall);
ucl_copy(this->tag, host_tag_write, false);
}
UCL_H_Vec<int> host_write;
host_write.view(sametag, max_same, *(this->ucl_device));
this->atom_sametag.resize_ib(max_same);
ucl_copy(this->atom_sametag, host_write, false);
this->map_array.resize_ib(map_size);
host_write.view(map_array, map_size, *(this->ucl_device));
ucl_copy(this->map_array, host_write, false);
}
}
// ---------------------------------------------------------------------------
// Copy nbor list from host if necessary and then calculate forces, virials,..
// ---------------------------------------------------------------------------
template <class numtyp, class acctyp>
void LJTIP4PLongT::compute(const int f_ago, const int inum_full,
const int nall, double **host_x, int *host_type,
int *ilist, int *numj, int **firstneigh,
const bool eflag_in, const bool vflag_in,
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) {
this->acc_timers();
int eflag, vflag;
if (eflag_in) eflag=2;
else eflag=0;
if (vflag_in) vflag=2;
else vflag=0;
this->set_kernel(eflag,vflag);
if (inum_full==0) {
host_start=0;
// Make sure textures are correct if realloc by a different hybrid style
this->resize_atom(0,nall,success);
this->zero_timers();
return;
}
int ago=this->hd_balancer.ago_first(f_ago);
int inum=this->hd_balancer.balance(ago,inum_full,cpu_time);
this->ans->inum(inum);
host_start=inum;
if (ago==0) {
this->reset_nbors(nall, inum, ilist, numj, firstneigh, success);
if (!success)
return;
}
this->atom->cast_x_data(host_x,host_type);
this->atom->cast_q_data(host_q);
this->hd_balancer.start_timer();
this->atom->add_x_data(host_x,host_type);
this->atom->add_q_data();
this->device->precompute(f_ago,nlocal,nall,host_x,host_type,success,host_q,
boxlo, prd);
t_ago = ago;
loop(eflag,vflag);
this->ans->copy_answers(eflag_in,vflag_in,eatom,vatom,ilist,inum);
this->device->add_ans_object(this->ans);
this->hd_balancer.stop_timer();
}
// ---------------------------------------------------------------------------
// Reneighbor on GPU if necessary and then compute forces, virials, energies
// ---------------------------------------------------------------------------
template <class numtyp, class acctyp>
int** LJTIP4PLongT::compute(const int ago, const int inum_full,
const int nall, double **host_x, int *host_type,
double *sublo, double *subhi, tagint *tag,
int *map_array, int map_size, int *sametag,
int max_same, int **nspecial, tagint **special,
const bool eflag_in, const bool vflag_in,
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) {
this->acc_timers();
int eflag, vflag;
if (eflag_in) eflag=2;
else eflag=0;
if (vflag_in) vflag=2;
else vflag=0;
this->set_kernel(eflag,vflag);
if (inum_full==0) {
host_start=0;
// Make sure textures are correct if realloc by a different hybrid style
this->resize_atom(0,nall,success);
this->zero_timers();
return nullptr;
}
this->hd_balancer.balance(cpu_time);
int inum=this->hd_balancer.get_gpu_count(ago,inum_full);
this->ans->inum(inum);
host_start=inum;
// Build neighbor list on GPU if necessary
if (ago==0) {
this->build_nbor_list(inum, inum_full-inum, nall, host_x, host_type,
sublo, subhi, tag, nspecial, special, success);
if (!success)
return nullptr;
this->atom->cast_q_data(host_q);
this->hd_balancer.start_timer();
} else {
this->atom->cast_x_data(host_x,host_type);
this->atom->cast_q_data(host_q);
this->hd_balancer.start_timer();
this->atom->add_x_data(host_x,host_type);
}
this->atom->add_q_data();
*ilist=this->nbor->host_ilist.begin();
*jnum=this->nbor->host_acc.begin();
copy_relations_data(nall, tag, map_array, map_size, sametag, max_same, ago);
this->device->precompute(ago,inum_full,nall,host_x,host_type,success,host_q,
boxlo, prd);
t_ago = ago;
loop(eflag,vflag);
this->ans->copy_answers(eflag_in,vflag_in,eatom,vatom,inum);
this->device->add_ans_object(this->ans);
this->hd_balancer.stop_timer();
return this->nbor->host_jlist.begin()-host_start;
}
template class LJ_TIP4PLong<PRECISION,ACC_PRECISION>;
}
|