Added Kokkos accelerated SLLOD thermostat (nvt/sllod/kk).

This commit is contained in:
Emily Kahl
2021-07-29 09:54:15 +10:00
parent e7ba4179a7
commit 8945d81be3
2 changed files with 268 additions and 0 deletions

View File

@ -0,0 +1,171 @@
// clang-format off
/* ----------------------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
https://www.lammps.org/
Steve Plimpton, sjplimp@sandia.gov, Sandia National Laboratories
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: Emily Kahl (UQ, e.kahl@uq.edu.au)
------------------------------------------------------------------------- */
#include "fix_nvt_sllod_kokkos.h"
#include "atom.h"
#include "compute.h"
#include "domain.h"
#include "error.h"
#include "fix.h"
#include "fix_deform_kokkos.h"
#include "group.h"
#include "math_extra.h"
#include "modify.h"
#include "atom_kokkos.h"
#include "atom.h"
#include "atom_masks.h"
#include "kokkos_few.h"
#include "memory_kokkos.h"
#include <cstring>
using namespace LAMMPS_NS;
using namespace FixConst;
/* ---------------------------------------------------------------------- */
template<class DeviceType>
FixNVTSllodKokkos<DeviceType>::FixNVTSllodKokkos(LAMMPS *lmp, int narg, char **arg) :
FixNHKokkos<DeviceType>(lmp, narg, arg)
{
/************************ EVK *********************/
atomKK = (AtomKokkos *) this->atom;
this->kokkosable = 1;
this->domainKK = (DomainKokkos *) this->domain;
if (!this->tstat_flag)
this->error->all(FLERR,"Temperature control must be used with fix nvt/kk");
if (this->pstat_flag)
this->error->all(FLERR,"Pressure control can not be used with fix nvt/kk");
if (this->mtchain_default_flag) this->mtchain = 1;
this->id_temp = utils::strdup(std::string(this->id)+"_temp");
this->modify->add_compute(fmt::format("{} all temp/deform/kk",this->id_temp));
this->tcomputeflag = 1;
}
/* ---------------------------------------------------------------------- */
template<class DeviceType>
FixNVTSllodKokkos<DeviceType>::~FixNVTSllodKokkos()
{
}
/* ---------------------------------------------------------------------- */
template<class DeviceType>
void FixNVTSllodKokkos<DeviceType>::init()
{
FixNHKokkos<DeviceType>::init();
vdelu = typename ArrayTypes<DeviceType>::t_v_array("nvt/sllod/kk:vdelu", atomKK->nlocal);
if (!this->temperature->tempbias)
this->error->all(FLERR,"Temperature for fix nvt/sllod does not have a bias");
nondeformbias = 0;
if (strncmp(this->temperature->style,"temp/deform", 11) != 0) nondeformbias = 1;
// check fix deform remap settings
int i;
for (i = 0; i < this->modify->nfix; i++)
if (strncmp(this->modify->fix[i]->style,"deform",6) == 0) {
if (((FixDeform *) this->modify->fix[i])->remapflag != Domain::V_REMAP)
this->error->all(FLERR,"Using fix nvt/sllod with inconsistent fix deform "
"remap option");
break;
}
if (i == this->modify->nfix)
this->error->all(FLERR,"Using fix nvt/sllod with no fix deform defined");
}
/* ----------------------------------------------------------------------
perform half-step scaling of velocities
-----------------------------------------------------------------------*/
template<class DeviceType>
void FixNVTSllodKokkos<DeviceType>::nh_v_temp()
{
// remove and restore bias = streaming velocity = Hrate*lamda + Hratelo
// thermostat thermal velocity only
// vdelu = SLLOD correction = Hrate*Hinv*vthermal
// for non temp/deform BIAS:
// calculate temperature since some computes require temp
// computed on current nlocal atoms to remove bias
if (nondeformbias){
atomKK->sync(this->temperature->execution_space,this->temperature->datamask_read);
this->temperature->compute_scalar();
atomKK->modified(this->temperature->execution_space,this->temperature->datamask_modify);
}
v = atomKK->k_v.view<DeviceType>();
mask = atomKK->k_mask.view<DeviceType>();
int nlocal = atomKK->nlocal;
if (this->igroup == atomKK->firstgroup) nlocal = atomKK->nfirst;
double h_two[6];
MathExtra::multiply_shape_shape(this->domain->h_rate,this->domain->h_inv,h_two);
d_h_two = Few<double, 6>(h_two);
this->copymode = 1;
Kokkos::parallel_for(Kokkos::RangePolicy<DeviceType, TagFixNVTSllod_temp1>(0,nlocal),*this);
this->copymode = 0;
this->temperature->remove_bias_all();
this->copymode = 1;
Kokkos::parallel_for(Kokkos::RangePolicy<DeviceType, TagFixNVTSllod_temp2>(0,nlocal),*this);
this->copymode = 0;
this->temperature->restore_bias_all();
}
template<class DeviceType>
KOKKOS_INLINE_FUNCTION
void FixNVTSllodKokkos<DeviceType>::operator()(TagFixNVTSllod_temp1, const int &i) const {
if (mask[i] & this->groupbit) {
vdelu(i,0) = d_h_two[0]*v(i,0) + d_h_two[5]*v(i,1) + d_h_two[4]*v(i,2);
vdelu(i,1) = d_h_two[1]*v(i,1) + d_h_two[3]*v(i,2);
vdelu(i,2) = d_h_two[2]*v(i,2);
}
}
template<class DeviceType>
KOKKOS_INLINE_FUNCTION
void FixNVTSllodKokkos<DeviceType>::operator()(TagFixNVTSllod_temp2, const int &i) const {
if (mask[i] & this->groupbit) {
v(i,0) = v(i,0)*this->factor_eta - this->dthalf*vdelu(i,0);
v(i,1) = v(i,1)*this->factor_eta - this->dthalf*vdelu(i,1);
v(i,2) = v(i,2)*this->factor_eta - this->dthalf*vdelu(i,2);
}
}
namespace LAMMPS_NS {
template class FixNVTSllodKokkos<LMPDeviceType>;
#ifdef LMP_KOKKOS_GPU
template class FixNVTSllodKokkos<LMPHostType>;
#endif
}

View File

@ -0,0 +1,97 @@
/* -*- c++ -*- ----------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
www.cs.sandia.gov/~sjplimp/lammps.html
Steve Plimpton, sjplimp@sandia.gov, Sandia National Laboratories
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
// clang-format off
FixStyle(nvt/sllod/kk,FixNVTSllodKokkos<LMPDeviceType>);
FixStyle(nvt/sllod/kk/device,FixNVTSllodKokkos<LMPDeviceType>);
FixStyle(nvt/sllod/kk/host,FixNVTSllodKokkos<LMPHostType>);
// clang-format on
#else
#ifndef LMP_FIX_NVT_SLLOD_KOKKOS_H
#define LMP_FIX_NVT_SLLOD_KOKKOS_H
#include "fix_nh_kokkos.h"
//#include "fix_nvt_sllod.h"
#include "kokkos_type.h"
#include "kokkos_few.h"
namespace LAMMPS_NS {
struct TagFixNVTSllod_temp1{};
struct TagFixNVTSllod_temp2{};
template<class DeviceType>
class FixNVTSllodKokkos : public FixNHKokkos<DeviceType> {
public:
FixNVTSllodKokkos(class LAMMPS *, int, char **);
~FixNVTSllodKokkos();
void init();
KOKKOS_INLINE_FUNCTION
void operator()(TagFixNVTSllod_temp1, const int& i) const;
KOKKOS_INLINE_FUNCTION
void operator()(TagFixNVTSllod_temp2, const int& i) const;
private:
int nondeformbias;
void nh_v_temp();
protected:
typename ArrayTypes<DeviceType>::t_x_array x;
typename ArrayTypes<DeviceType>::t_v_array v;
typename ArrayTypes<DeviceType>::t_v_array vdelu;
typename ArrayTypes<DeviceType>::t_f_array_const f;
typename ArrayTypes<DeviceType>::t_float_1d rmass;
typename ArrayTypes<DeviceType>::t_float_1d mass;
typename ArrayTypes<DeviceType>::t_int_1d type;
typename ArrayTypes<DeviceType>::t_int_1d mask;
Few<double, 6> d_h_two;
class DomainKokkos *domainKK;
class AtomKokkos *atomKK;
};
} // namespace LAMMPS_NS
#endif
#endif
/* ERROR/WARNING messages:
E: Temperature control must be used with fix nvt/sllod
Self-explanatory.
E: Pressure control can not be used with fix nvt/sllod
Self-explanatory.
E: Temperature for fix nvt/sllod does not have a bias
The specified compute must compute temperature with a bias.
E: Using fix nvt/sllod with inconsistent fix deform remap option
Fix nvt/sllod requires that deforming atoms have a velocity profile
provided by "remap v" as a fix deform option.
E: Using fix nvt/sllod with no fix deform defined
Self-explanatory.
*/