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: Tim Fuller (Sandia)
17 ------------------------------------------------------------------------- */
18
19 #include "mpi.h"
20 #include "string.h"
21 #include "compute_temp_kokkos.h"
22 #include "update.h"
23 #include "particle_kokkos.h"
24 #include "domain.h"
25 #include "error.h"
26 #include "kokkos.h"
27 #include "sparta_masks.h"
28
29 using namespace SPARTA_NS;
30
31 /* ---------------------------------------------------------------------- */
32
ComputeTempKokkos(SPARTA * sparta,int narg,char ** arg)33 ComputeTempKokkos::ComputeTempKokkos(SPARTA *sparta, int narg, char **arg) :
34 ComputeTemp(sparta, narg, arg)
35 {
36 if (narg != 2) error->all(FLERR,"Illegal compute temp command");
37 scalar_flag = 1;
38 kokkos_flag = 1;
39 }
40
41 /* ---------------------------------------------------------------------- */
42
compute_scalar()43 double ComputeTempKokkos::compute_scalar()
44 {
45 if (sparta->kokkos->prewrap) {
46 return ComputeTemp::compute_scalar();
47 } else {
48 return compute_scalar_kokkos();
49 }
50 }
51
52 KOKKOS_INLINE_FUNCTION
operator ()(const int & i,double & lsum) const53 void ComputeTempKokkos::operator()(const int& i, double& lsum) const {
54 double* v = d_particles[i].v;
55 const int ispecies = d_particles[i].ispecies;
56 const double mass = d_species[ispecies].mass;
57 lsum += (v[0]*v[0] + v[1]*v[1] + v[2]*v[2]) * mass;
58 }
59
compute_scalar_kokkos()60 double ComputeTempKokkos::compute_scalar_kokkos()
61 {
62 copymode = 1;
63 invoked_scalar = update->ntimestep;
64 ParticleKokkos* particle_kk = (ParticleKokkos*) particle;
65 particle_kk->sync(Device, PARTICLE_MASK|SPECIES_MASK);
66 d_particles = particle_kk->k_particles.d_view;
67 d_species = particle_kk->k_species.d_view;
68 int nlocal = particle->nlocal;
69
70 double t = 0.0;
71 auto range_policy = Kokkos::RangePolicy<DeviceType>(0, nlocal);
72 Kokkos::parallel_reduce(range_policy, *this, t);
73 DeviceType().fence();
74 copymode = 0;
75
76 d_particles = t_particle_1d(); // destroy reference to reduce memory use
77
78 MPI_Allreduce(&t,&scalar,1,MPI_DOUBLE,MPI_SUM,world);
79
80 bigint n = particle->nlocal;
81 MPI_Allreduce(&n,&particle->nglobal,1,MPI_SPARTA_BIGINT,MPI_SUM,world);
82 if (particle->nglobal == 0) return 0.0;
83
84 // normalize with 3 instead of dim since even 2d has 3 velocity components
85
86 double factor = update->mvv2e / (3.0 * particle->nglobal * update->boltz);
87 scalar *= factor;
88 return scalar;
89 }
90