diff options
author | 2022-06-01 17:34:11 -0700 | |
---|---|---|
committer | 2022-06-01 17:34:11 -0700 | |
commit | 5e65abc3950b8c08422a15073d17e30bb8c0e80b (patch) | |
tree | 0f0f6afaa480dc87389c5df6308015fa523c4e67 /Source/Python/WarpXWrappers.cpp | |
parent | 40860cf27d9d2b99f693da384aaa63e6642192a4 (diff) | |
download | WarpX-5e65abc3950b8c08422a15073d17e30bb8c0e80b.tar.gz WarpX-5e65abc3950b8c08422a15073d17e30bb8c0e80b.tar.zst WarpX-5e65abc3950b8c08422a15073d17e30bb8c0e80b.zip |
Clean-up in python wrappers to access particle data (#2531)
* use numLocalTilesAtLevel instead of an empty particle iterator loop to determine the number of non-empty tiles
* explicitly zero out the `*particles_per_tile` array
Diffstat (limited to 'Source/Python/WarpXWrappers.cpp')
-rw-r--r-- | Source/Python/WarpXWrappers.cpp | 45 |
1 files changed, 15 insertions, 30 deletions
diff --git a/Source/Python/WarpXWrappers.cpp b/Source/Python/WarpXWrappers.cpp index cb88e3a89..399597f60 100644 --- a/Source/Python/WarpXWrappers.cpp +++ b/Source/Python/WarpXWrappers.cpp @@ -446,15 +446,12 @@ namespace const std::string species_name(char_species_name); auto & myspc = mypc.GetParticleContainerFromName(species_name); - int i = 0; - for (WarpXParIter pti(myspc, lev); pti.isValid(); ++pti, ++i) {} - - // *num_tiles = myspc.numLocalTilesAtLevel(lev); - *num_tiles = i; + *num_tiles = myspc.numLocalTilesAtLevel(lev); *particles_per_tile = static_cast<int*>(malloc(*num_tiles*sizeof(int))); + memset(*particles_per_tile, 0, *num_tiles*sizeof(int)); auto data = static_cast<amrex::ParticleReal**>(malloc(*num_tiles*sizeof(typename WarpXParticleContainer::ParticleType*))); - i = 0; + int i = 0; for (WarpXParIter pti(myspc, lev); pti.isValid(); ++pti, ++i) { auto& aos = pti.GetArrayOfStructs(); data[i] = (amrex::ParticleReal*) aos.data(); @@ -473,15 +470,12 @@ namespace int comp = warpx_getParticleCompIndex(char_species_name, char_comp_name); - int i = 0; - for (WarpXParIter pti(myspc, lev); pti.isValid(); ++pti, ++i) {} - - // *num_tiles = myspc.numLocalTilesAtLevel(lev); - *num_tiles = i; + *num_tiles = myspc.numLocalTilesAtLevel(lev); *particles_per_tile = static_cast<int*>(malloc(*num_tiles*sizeof(int))); + memset(*particles_per_tile, 0, *num_tiles*sizeof(int)); auto data = static_cast<amrex::ParticleReal**>(malloc(*num_tiles*sizeof(amrex::ParticleReal*))); - i = 0; + int i = 0; for (WarpXParIter pti(myspc, lev); pti.isValid(); ++pti, ++i) { auto& soa = pti.GetStructOfArrays(); data[i] = (amrex::ParticleReal*) soa.GetRealData(comp).dataPtr(); @@ -541,15 +535,12 @@ namespace const int comp = particle_buffer.NumIntComps() - 1; - int i = 0; - for (amrex::ParIter<0,0,PIdx::nattribs, 0, amrex::PinnedArenaAllocator> pti(particle_buffer, lev); pti.isValid(); ++pti, ++i) {} - - // *num_tiles = myspc.numLocalTilesAtLevel(lev); - *num_tiles = i; + *num_tiles = particle_buffer.numLocalTilesAtLevel(lev); *particles_per_tile = static_cast<int*>(malloc(*num_tiles*sizeof(int))); + memset(*particles_per_tile, 0, *num_tiles*sizeof(int)); auto data = static_cast<int**>(malloc(*num_tiles*sizeof(int*))); - i = 0; + int i = 0; for (amrex::ParIter<0,0,PIdx::nattribs, 0, amrex::PinnedArenaAllocator> pti(particle_buffer, lev); pti.isValid(); ++pti, ++i) { auto& soa = pti.GetStructOfArrays(); data[i] = (int*) soa.GetIntData(comp).dataPtr(); @@ -568,15 +559,12 @@ namespace const int comp = warpx_getParticleCompIndex(species_name, comp_name); - int i = 0; - for (amrex::ParIter<0,0,PIdx::nattribs, 0, amrex::PinnedArenaAllocator> pti(particle_buffer, lev); pti.isValid(); ++pti, ++i) {} - - // *num_tiles = myspc.numLocalTilesAtLevel(lev); - *num_tiles = i; + *num_tiles = particle_buffer.numLocalTilesAtLevel(lev); *particles_per_tile = static_cast<int*>(malloc(*num_tiles*sizeof(int))); + memset(*particles_per_tile, 0, *num_tiles*sizeof(int)); auto data = static_cast<amrex::ParticleReal**>(malloc(*num_tiles*sizeof(amrex::ParticleReal*))); - i = 0; + int i = 0; for (amrex::ParIter<0,0,PIdx::nattribs, 0, amrex::PinnedArenaAllocator> pti(particle_buffer, lev); pti.isValid(); ++pti, ++i) { auto& soa = pti.GetStructOfArrays(); data[i] = (amrex::ParticleReal*) soa.GetRealData(comp).dataPtr(); @@ -593,15 +581,12 @@ namespace auto& particle_buffers = WarpX::GetInstance().GetParticleBoundaryBuffer(); auto& particle_buffer = particle_buffers.getParticleBuffer(species_name, boundary); - int i = 0; - for (amrex::ParIter<0,0,PIdx::nattribs, 0, amrex::PinnedArenaAllocator> pti(particle_buffer, lev); pti.isValid(); ++pti, ++i) {} - - // *num_tiles = myspc.numLocalTilesAtLevel(lev); - *num_tiles = i; + *num_tiles = particle_buffer.numLocalTilesAtLevel(lev); *particles_per_tile = static_cast<int*>(malloc(*num_tiles*sizeof(int))); + memset(*particles_per_tile, 0, *num_tiles*sizeof(int)); auto data = static_cast<amrex::ParticleReal**>(malloc(*num_tiles*sizeof(typename WarpXParticleContainer::ParticleType*))); - i = 0; + int i = 0; for (amrex::ParIter<0,0,PIdx::nattribs, 0, amrex::PinnedArenaAllocator> pti(particle_buffer, lev); pti.isValid(); ++pti, ++i) { auto& aos = pti.GetArrayOfStructs(); data[i] = (amrex::ParticleReal*) aos.data(); |