Files
lammps/lib/gpu/lal_gran_hooke.cpp

393 lines
13 KiB
C++

/***************************************************************************
* gran_hooke.cpp
* -------------------
* Trung Dac Nguyen (ORNL)
*
* Class for acceleration of the gran/hooke pair style.
*
* __________________________________________________________________________
* This file is part of the LAMMPS Accelerator Library (LAMMPS_AL)
* __________________________________________________________________________
*
* begin :
* email : nguyentd@ornl.gov
***************************************************************************/
#if defined(USE_OPENCL)
#include "gran_hooke_cl.h"
#elif defined(USE_CUDART)
const char *gran_hooke=0;
#else
#include "gran_hooke_cubin.h"
#endif
#include "lal_gran_hooke.h"
#include <cassert>
using namespace LAMMPS_AL;
#define GranHookeT GranHooke<numtyp, acctyp>
extern Device<PRECISION,ACC_PRECISION> device;
template <class numtyp, class acctyp>
GranHookeT::GranHooke() : BaseSphere<numtyp,acctyp>(),
_allocated(false) {
}
template <class numtyp, class acctyp>
GranHookeT::~GranHooke() {
clear();
}
template <class numtyp, class acctyp>
int GranHookeT::init(const int ntypes,
double host_k_n, double host_k_t,
double host_gamman, double host_gammat,
double host_xmu, double host_dt, const int host_dampflag,
double *host_special_lj, int *host_mask,
const int nlocal, const int nall, const int max_nbors,
const int maxspecial, const double cell_size,
const double gpu_split, FILE *screen) {
int onetype=0;
int success;
success = this->init_atomic(nlocal, nall, max_nbors, maxspecial,
cell_size, gpu_split, screen,
gran_hooke, "k_gran_hooke");
if (success != 0)
return success;
_shared_view = this->ucl_device->shared_memory() && sizeof(numtyp)==sizeof(double);
// allocate rad
int ef_nall=nall;
if (ef_nall==0)
ef_nall=50000;
_max_rad_size=static_cast<int>(static_cast<double>(ef_nall)*1.10);
if (!_shared_view)
{
c_rad.alloc(_max_rad_size,*(this->ucl_device),UCL_WRITE_ONLY,UCL_READ_ONLY);
c_rmass.alloc(_max_rad_size,*(this->ucl_device),UCL_WRITE_ONLY,UCL_READ_ONLY);
c_omega.alloc(_max_rad_size*4,*(this->ucl_device),UCL_WRITE_ONLY,UCL_READ_ONLY);
c_vel.alloc(_max_rad_size*4,*(this->ucl_device),UCL_WRITE_ONLY,UCL_READ_ONLY);
}
rad_tex.get_texture(*(this->pair_program),"rad_tex");
rad_tex.bind_float(c_rad,1);
rmass_tex.get_texture(*(this->pair_program),"rmass_tex");
rmass_tex.bind_float(c_rmass,1);
omega_tex.get_texture(*(this->pair_program),"omega_tex");
omega_tex.bind_float(c_omega,4);
vel_tex.get_texture(*(this->pair_program),"vel_tex");
vel_tex.bind_float(c_vel,4);
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) {//TODO
lj_types=max_shared_types;
shared_types=true;
}
//shared_types=false;//TEST
_k_n=host_k_n;
_k_t=host_k_t;
_gamman=host_gamman;
_gammat=host_gammat;
_xmu=host_xmu;
_dampflag=host_dampflag;
_dt=host_dt;
UCL_H_Vec<int> h_mask;
_mask.alloc(_max_rad_size,*(this->ucl_device),UCL_READ_ONLY);
h_mask.view(host_mask,_max_rad_size,*(this->ucl_device));
ucl_copy(_mask,h_mask,false);
UCL_H_Vec<double> dview;
sp_lj.alloc(4,*(this->ucl_device),UCL_READ_ONLY);
dview.view(host_special_lj,4,*(this->ucl_device));
ucl_copy(sp_lj,dview,false);
_allocated=true;
this->_max_bytes=sp_lj.row_bytes() + _mask.row_bytes();
return 0;
}
template <class numtyp, class acctyp>
void GranHookeT::clear() {
if (!_allocated)
return;
_allocated = false;
sp_lj.clear();
_mask.clear();
c_rad.clear();
c_omega.clear();
c_rmass.clear();
c_vel.clear();
this->clear_atomic();
}
template <class numtyp, class acctyp>
double GranHookeT::host_memory_usage() const {
return this->host_memory_usage_atomic() +
sizeof(GranHooke<numtyp,acctyp>);
}
template <class numtyp, class acctyp>
int GranHookeT::bytes_per_atom(const int max_nbors) const {
return this->bytes_per_atom_atomic(max_nbors) +
sizeof(numtyp4)*2 + sizeof(numtyp2) + sizeof(int) +
sizeof(numtyp)*3 + sizeof(numtyp)*3*max_nbors;
}
// ---------------------------------------------------------------------------
// Calculate energies, forces, and torques
// ---------------------------------------------------------------------------
template <class numtyp, class acctyp>
int GranHookeT::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 GX = static_cast<int>(ceil(static_cast<double>(this->ans->inum()) / BX));
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();
if (shared_types) {
this->k_pair_fast.set_size(GX, BX);
this->k_pair_fast.run(&this->atom->x, &c_rad, &c_vel,
&c_omega, &sp_lj, &c_rmass,
&this->nbor->dev_nbor, &this->_nbor_data->begin(),
&_mask, &this->ans->force,
&this->ans->engv, &_eflag, &_vflag,
&ainum, &nbor_pitch, &this->_threads_per_atom,
&_k_n, &_k_t, &_gamman, &_gammat, &_xmu,
&_dampflag, &freeze_group_bit,
&_limit_damping, &_dt);
} else {
this->k_pair.set_size(GX, BX);
this->k_pair.run(&this->atom->x, &c_rad, &c_vel,
&c_omega, &sp_lj, &c_rmass,
&this->nbor->dev_nbor, &this->_nbor_data->begin(),
&_mask, &this->ans->force,
&this->ans->engv, &_eflag, &_vflag,
&ainum, &nbor_pitch, &this->_threads_per_atom,
&_k_n, &_k_t, &_gamman, &_gammat, &_xmu,
&_dampflag, &freeze_group_bit,
&_limit_damping, &_dt);
}
this->time_pair.stop();
return GX;
}
template <class numtyp, class acctyp>
void GranHookeT::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, tagint *tag,
double **host_v, double *host_rad,
double **host_omega, double *host_rmass, int limit_damping,
const int nlocal, double *boxlo, double *prd) {
this->acc_timers();
int eflag, vflag;
if (eatom) eflag=2;
else if (eflag_in) eflag=1;
else eflag=0;
if (vatom) vflag=2;
else if (vflag_in) vflag=1;
else vflag=0;
#ifdef LAL_NO_BLOCK_REDUCE
if (eflag) eflag=2;
if (vflag) vflag=2;
#endif
this->set_kernel(eflag,vflag);
_limit_damping = limit_damping;
// ------------------- Resize data array --------------------------
if (nall>_max_rad_size) {
_max_rad_size=static_cast<int>(static_cast<double>(nall)*1.10);
if (!_shared_view) {
c_rad.resize(_max_rad_size);
rad_tex.bind_float(c_rad,1);
c_omega.resize(_max_rad_size*4);
omega_tex.bind_float(c_omega,4);
c_vel.resize(_max_rad_size*4);
vel_tex.bind_float(c_vel,4);
c_rmass.resize(_max_rad_size);
rmass_tex.bind_float(c_rmass,1);
}
}
// ----------------------------------------------------------------
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;
}
// load balance, returning the atom count on the device (inum)
int ago=this->hd_balancer.ago_first(f_ago);
//this->hd_balancer.balance(cpu_time);
//int inum=this->hd_balancer.get_gpu_count(ago,inum_full);
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->cast_vel_data(host_v);
this->cast_rad_data(host_rad);
this->cast_omega_data(host_omega);
this->cast_rmass_data(host_rmass);
this->hd_balancer.start_timer();
this->atom->add_x_data(host_x,host_type);
this->add_rad_data();
this->add_vel_data();
this->add_rmass_data();
this->add_omega_data();
const int red_blocks=this->loop(eflag,vflag);
this->ans->copy_answers(eflag_in,vflag_in,eatom,vatom,ilist,red_blocks);
this->device->add_ans_object(this->ans);
this->hd_balancer.stop_timer();
}
template <class numtyp, class acctyp>
int** GranHookeT::compute(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_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_v, double *host_rad,
double **host_omega, double *host_rmass, int limit_damping,
double *boxlo, double *prd) {
this->acc_timers();
int eflag, vflag;
if (eatom) eflag=2;
else if (eflag_in) eflag=1;
else eflag=0;
if (vatom) vflag=2;
else if (vflag_in) vflag=1;
else vflag=0;
#ifdef LAL_NO_BLOCK_REDUCE
if (eflag) eflag=2;
if (vflag) vflag=2;
#endif
this->set_kernel(eflag,vflag);
_limit_damping= limit_damping;
// ------------------- Resize data array --------------------------
if (nall>_max_rad_size) {
_max_rad_size=static_cast<int>(static_cast<double>(nall)*1.10);
if (!_shared_view) {
c_rad.resize(_max_rad_size);
rad_tex.bind_float(c_rad,1);
c_omega.resize(_max_rad_size * 4);
omega_tex.bind_float(c_omega,4);
c_vel.resize(_max_rad_size * 4);
vel_tex.bind_float(c_vel,4);
c_rmass.resize(_max_rad_size);
rmass_tex.bind_float(c_rmass,1);
}
}
// -----------------------------------------------------------------
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;
}
// ----------------------------------------------------------------
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;
}
// load balance, returning the atom count on the device (inum)
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->cast_rad_data(host_rad);
this->cast_vel_data(host_v);
this->cast_omega_data(host_omega);
this->cast_rmass_data(host_rmass);
this->hd_balancer.start_timer();
} else {
this->atom->cast_x_data(host_x,host_type);
this->cast_vel_data(host_v);
this->cast_rad_data(host_rad);
this->cast_omega_data(host_omega);
this->cast_rmass_data(host_rmass);
this->hd_balancer.start_timer();
this->atom->add_x_data(host_x,host_type);
}
this->add_rad_data();
this->add_omega_data();
this->add_vel_data();
this->add_rmass_data();
*ilist=this->nbor->host_ilist.begin();
*jnum=this->nbor->host_acc.begin();
const int red_blocks=this->loop(eflag,vflag);
this->ans->copy_answers(eflag_in,vflag_in,eatom,vatom,red_blocks);
this->device->add_ans_object(this->ans);
this->hd_balancer.stop_timer();
return this->nbor->host_jlist.begin()-host_start;
}
// Instantiate the class
template class GranHooke<PRECISION,ACC_PRECISION>;