Adding fix_wall_lj93_kokkos

This commit is contained in:
Stan Moore
2017-03-02 14:02:49 -07:00
parent 8210b25fb8
commit 3820c5881d
5 changed files with 193 additions and 2 deletions

View File

@ -105,6 +105,8 @@ action fix_setforce_kokkos.cpp
action fix_setforce_kokkos.h
action fix_momentum_kokkos.cpp
action fix_momentum_kokkos.h
action fix_wall_lj93_kokkos.cpp
action fix_wall_lj93_kokkos.h
action fix_wall_reflect_kokkos.cpp
action fix_wall_reflect_kokkos.h
action fix_dpd_energy_kokkos.cpp fix_dpd_energy.cpp

View File

@ -0,0 +1,104 @@
/* ----------------------------------------------------------------------
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.
------------------------------------------------------------------------- */
#include <math.h>
#include "fix_wall_lj93_kokkos.h"
#include "atom_kokkos.h"
#include "error.h"
#include "atom_masks.h"
using namespace LAMMPS_NS;
using namespace FixConst;
/* ---------------------------------------------------------------------- */
template <class DeviceType>
FixWallLJ93Kokkos<DeviceType>::FixWallLJ93Kokkos(LAMMPS *lmp, int narg, char **arg) :
FixWallLJ93(lmp, narg, arg)
{
kokkosable = 1;
atomKK = (AtomKokkos *) atom;
execution_space = ExecutionSpaceFromDevice<DeviceType>::space;
datamask_read = EMPTY_MASK;
datamask_modify = EMPTY_MASK;
}
/* ----------------------------------------------------------------------
interaction of all particles in group with a wall
m = index of wall coeffs
which = xlo,xhi,ylo,yhi,zlo,zhi
error if any particle is on or behind wall
------------------------------------------------------------------------- */
template <class DeviceType>
void FixWallLJ93Kokkos<DeviceType>::wall_particle(int m_in, int which, double coord_in)
{
m = m_in;
coord = coord_in;
atomKK->sync(execution_space, X_MASK|F_MASK|MASK_MASK);
x = atomKK->k_x.view<DeviceType>();
f = atomKK->k_f.view<DeviceType>();
mask = atomKK->k_mask.view<DeviceType>();
DAT::tdual_int_scalar k_oneflag = DAT::tdual_int_scalar("fix:oneflag");
d_oneflag = k_oneflag.view<DeviceType>();
int nlocal = atom->nlocal;
dim = which / 2;
side = which % 2;
if (side == 0) side = -1;
copymode = 1;
FixWallLJ93KokkosFunctor<DeviceType> wp_functor(this);
Kokkos::parallel_reduce(nlocal,wp_functor,ewall);
DeviceType::fence();
copymode = 0;
atomKK->modified(execution_space, F_MASK);
k_oneflag.template modify<DeviceType>();
k_oneflag.template sync<LMPHostType>();
if (k_oneflag.h_view()) error->one(FLERR,"Particle on or inside fix wall surface");
}
template <class DeviceType>
KOKKOS_INLINE_FUNCTION
void FixWallLJ93Kokkos<DeviceType>::wall_particle_item(int i, value_type ewall) const {
if (mask(i) & groupbit) {
double delta;
if (side < 0) delta = x(i,dim) - coord;
else delta = coord - x(i,dim);
if (delta >= cutoff[m]) return;
if (delta <= 0.0) {
d_oneflag() = 1;
return;
}
double rinv = 1.0/delta;
double r2inv = rinv*rinv;
double r4inv = r2inv*r2inv;
double r10inv = r4inv*r4inv*r2inv;
double fwall = side * (coeff1[m]*r10inv - coeff2[m]*r4inv);
f(i,dim) -= fwall;
ewall[0] += coeff3[m]*r4inv*r4inv*rinv -
coeff4[m]*r2inv*rinv - offset[m];
ewall[m+1] += fwall;
}
}
namespace LAMMPS_NS {
template class FixWallLJ93Kokkos<LMPDeviceType>;
#ifdef KOKKOS_HAVE_CUDA
template class FixWallLJ93Kokkos<LMPHostType>;
#endif
}

View File

@ -0,0 +1,83 @@
/* -*- c++ -*- ----------------------------------------------------------
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.
------------------------------------------------------------------------- */
#ifdef FIX_CLASS
FixStyle(wall/lj93/kk,FixWallLJ93Kokkos<LMPDeviceType>)
FixStyle(wall/lj93/kk/device,FixWallLJ93Kokkos<LMPDeviceType>)
FixStyle(wall/lj93/kk/host,FixWallLJ93Kokkos<LMPHostType>)
#else
#ifndef LMP_FIX_WALL_LJ93_KOKKOS_H
#define LMP_FIX_WALL_LJ93_KOKKOS_H
#include "fix_wall_lj93.h"
#include "kokkos_type.h"
namespace LAMMPS_NS {
template <class DeviceType>
class FixWallLJ93Kokkos : public FixWallLJ93 {
public:
typedef DeviceType device_type;
typedef ArrayTypes<DeviceType> AT;
typedef double value_type[];
FixWallLJ93Kokkos(class LAMMPS *, int, char **);
void wall_particle(int, int, double);
int m;
KOKKOS_INLINE_FUNCTION
void wall_particle_item(int, value_type) const;
private:
int dim,side;
double coord;
typename AT::t_x_array x;
typename AT::t_f_array f;
typename AT::t_int_1d mask;
typename AT::t_int_scalar d_oneflag;
};
template <class DeviceType>
struct FixWallLJ93KokkosFunctor {
typedef DeviceType device_type ;
typedef double value_type[];
const int value_count;
FixWallLJ93Kokkos<DeviceType> c;
FixWallLJ93KokkosFunctor(FixWallLJ93Kokkos<DeviceType>* c_ptr):
c(*c_ptr),
value_count(c.m) {}
KOKKOS_INLINE_FUNCTION
void operator()(const int i, value_type ewall) const {
c.wall_particle_item(i,ewall);
}
};
}
#endif
#endif
/* ERROR/WARNING messages:
E: Particle on or inside fix wall surface
Particles must be "exterior" to the wall in order for energy/force to
be calculated.
*/

View File

@ -201,6 +201,8 @@ FixWall::FixWall(LAMMPS *lmp, int narg, char **arg) :
FixWall::~FixWall()
{
if (copymode) return;
for (int m = 0; m < nwall; m++) {
delete [] xstr[m];
delete [] estr[m];

View File

@ -28,9 +28,9 @@ class FixWallLJ93 : public FixWall {
public:
FixWallLJ93(class LAMMPS *, int, char **);
void precompute(int);
void wall_particle(int, int, double);
virtual void wall_particle(int, int, double);
private:
protected:
double coeff1[6],coeff2[6],coeff3[6],coeff4[6],offset[6];
};