1 /* ----------------------------------------------------------------------
2    SPARTA - Stochastic PArallel Rarefied-gas Time-accurate Analyzer
3    http://sparta.sandia.gov
4    Steve Plimpton, sjplimp@sandia.gov, Michael Gallis, magalli@sandia.gov
5    Sandia National Laboratories
6 
7    Copyright (2014) Sandia Corporation.  Under the terms of Contract
8    DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
9    certain rights in this software.  This software is distributed under
10    the GNU General Public License.
11 
12    See the README file in the top-level SPARTA directory.
13 ------------------------------------------------------------------------- */
14 
15 /* ----------------------------------------------------------------------
16    Contributing author: Alan Stagg (Sandia)
17 ------------------------------------------------------------------------- */
18 
19 #include "string.h"
20 #include "compute_ke_particle_kokkos.h"
21 #include "particle_kokkos.h"
22 #include "memory_kokkos.h"
23 #include "kokkos.h"
24 #include "update.h"
25 #include "sparta_masks.h"
26 #include "modify.h"
27 #include "comm.h"
28 #include "memory.h"
29 #include "error.h"
30 
31 using namespace SPARTA_NS;
32 
33 /* ---------------------------------------------------------------------- */
34 
ComputeKEParticleKokkos(SPARTA * sparta,int narg,char ** arg)35 ComputeKEParticleKokkos::ComputeKEParticleKokkos(SPARTA *sparta, int narg, char **arg) :
36   ComputeKEParticle(sparta, narg, arg)
37 {
38   kokkos_flag = 1;
39 }
40 
41 /* ---------------------------------------------------------------------- */
42 
~ComputeKEParticleKokkos()43 ComputeKEParticleKokkos::~ComputeKEParticleKokkos()
44 {
45   if (copymode) return;
46   memoryKK->destroy_kokkos(k_vector_particle,vector_particle);
47   vector_particle = NULL;
48 }
49 
50 /* ---------------------------------------------------------------------- */
51 
compute_per_particle()52 void ComputeKEParticleKokkos::compute_per_particle()
53 {
54   if (sparta->kokkos->prewrap) {
55     ComputeKEParticle::compute_per_particle();
56   }  else {
57     compute_per_particle_kokkos();
58   }
59 }
60 
61 /* ---------------------------------------------------------------------- */
62 
compute_per_particle_kokkos()63 void ComputeKEParticleKokkos::compute_per_particle_kokkos()
64 {
65   // grow ke array (d_vector) if necessary
66   if (particle->nlocal > nmax) {
67     memoryKK->destroy_kokkos(k_vector_particle,vector_particle);
68     nmax = particle->maxlocal;
69     memoryKK->create_kokkos(k_vector_particle,vector_particle,nmax,"ke/particle:vector_particle");
70     d_vector = k_vector_particle.d_view;
71   }
72 
73   ParticleKokkos* particle_kk = (ParticleKokkos*) particle;
74   particle_kk->sync(Device,PARTICLE_MASK|SPECIES_MASK);
75   d_particles = particle_kk->k_particles.d_view;
76   d_species = particle_kk->k_species.d_view;
77   int nlocal = particle->nlocal;
78 
79   // compute kinetic energy for each atom in group
80   copymode = 1;
81   Kokkos::parallel_for(Kokkos::RangePolicy<DeviceType>(0,nlocal),*this);
82   DeviceType().fence();
83   copymode = 0;
84 
85   d_particles = t_particle_1d(); // destroy reference to reduce memory use
86 }
87 
88 KOKKOS_INLINE_FUNCTION
operator ()(const int & i) const89 void ComputeKEParticleKokkos::operator()(const int &i) const {
90   const int ispecies = d_particles[i].ispecies;
91   const double mass = d_species[ispecies].mass;
92   const double mvv2e = update->mvv2e;
93   double *v = d_particles[i].v;
94   d_vector[i] = 0.5 * mvv2e * mass * (v[0]*v[0] + v[1]*v[1] + v[2]*v[2]);
95 }
96