1 #ifndef TRANSPOSE_PARTICLES_HPP
2 #define TRANSPOSE_PARTICLES_HPP
6 template<
class DataType>
10 if (VEC_LEN == 1)
return;
12 const std::size_t PTL_N_DBL = Particles::sizeof_ptl<DataType>()/8;
18 int team_size = VEC_LEN;
22 int league_size = n_vec;
24 typedef Kokkos::View<double[VEC_LEN*PTL_N_DBL],ExSpace::scratch_memory_space, Kokkos::MemoryTraits<Kokkos::Unmanaged>> PtlVec;
25 size_t shmem_size = PtlVec::shmem_size();
27 Kokkos::TeamPolicy<ExSpace> (league_size, team_size).set_scratch_size(
KOKKOS_TEAM_SCRATCH_OPT,Kokkos::PerTeam(shmem_size)),
28 KOKKOS_LAMBDA (Kokkos::TeamPolicy<ExSpace>::member_type team_member){
32 int ivec = team_member.league_rank();
33 int ith = team_member.team_rank();
36 for(
int idbl = ith; idbl<VEC_LEN*PTL_N_DBL; idbl+=team_size){
37 shmem_ptl_vec(idbl) = vptl[ivec].
data[idbl];
41 for(
int iptl = ith; iptl<VEC_LEN; iptl+=team_size){
42 for(
int iprop = 0; iprop<PTL_N_DBL; iprop++){
43 vptl[ivec].
data[iprop*VEC_LEN + iptl] = shmem_ptl_vec(iptl*PTL_N_DBL + iprop);
51 template<
class DataType>
54 if (VEC_LEN == 1)
return;
56 const std::size_t PTL_N_DBL = Particles::sizeof_ptl<DataType>()/8;
62 int team_size = VEC_LEN;
66 int league_size = n_vec;
68 typedef Kokkos::View<double[VEC_LEN*PTL_N_DBL],ExSpace::scratch_memory_space, Kokkos::MemoryTraits<Kokkos::Unmanaged>> PtlVec;
69 size_t shmem_size = PtlVec::shmem_size();
71 Kokkos::TeamPolicy<ExSpace> (league_size, team_size).set_scratch_size(
KOKKOS_TEAM_SCRATCH_OPT,Kokkos::PerTeam(shmem_size)),
72 KOKKOS_LAMBDA (Kokkos::TeamPolicy<ExSpace>::member_type team_member){
76 int ivec = team_member.league_rank();
77 int ith = team_member.team_rank();
80 for(
int iptl = ith; iptl<VEC_LEN; iptl+=team_size){
81 for(
int iprop = 0; iprop<PTL_N_DBL; iprop++){
82 shmem_ptl_vec(iptl*PTL_N_DBL + iprop) = vptl[ivec].
data[iprop*VEC_LEN + iptl];
87 for(
int idbl = ith; idbl<VEC_LEN*PTL_N_DBL; idbl+=team_size){
88 vptl[ivec].
data[idbl] = shmem_ptl_vec(idbl);
Definition: particles.hpp:248
void parallel_for(const std::string name, int n_ptl, Function func, Option option, HostAoSoA aosoa_h, DeviceAoSoA aosoa_d)
Definition: streamed_parallel_for.hpp:252
type(ptl_type), dimension(:), allocatable, target local_particles
Definition: resamp_mod.F90:150
#define KOKKOS_TEAM_SCRATCH_OPT
Definition: space_settings.hpp:77
Definition: particles.hpp:123
double data[T *VEC_LEN]
Definition: particles.hpp:124
void transpose_particles_from_AoS_to_AoSoA(Particles::Array< DataType, DeviceType > &local_particles, int ioffset, int n_vec)
Definition: transpose_particles.hpp:7
void transpose_particles_from_AoSoA_to_AoS(Particles::Array< DataType, DeviceType > &local_particles, int ioffset, int n_vec)
Definition: transpose_particles.hpp:52