Skip to content

Commit

Permalink
Update SortParticlesForDeposition for pure SoA (#3277)
Browse files Browse the repository at this point in the history
I also fixed a typo in the SoA particle `pos` routine that returned a
`RealVect`.

The proposed changes:
- [x] fix a bug or incorrect behavior in AMReX
- [x] add new capabilities to AMReX
- [ ] changes answers in the test suite to more than roundoff level
- [ ] are likely to significantly affect the results of downstream AMReX
users
- [ ] include documentation in the code and/or rst files, if appropriate
  • Loading branch information
atmyers authored Apr 26, 2023
1 parent ba69155 commit e7e2e1e
Show file tree
Hide file tree
Showing 3 changed files with 17 additions and 15 deletions.
8 changes: 3 additions & 5 deletions Src/Particle/AMReX_ParticleContainerI.H
Original file line number Diff line number Diff line change
Expand Up @@ -1203,16 +1203,14 @@ ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator>

for(MFIter mfi = MakeMFIter(lev); mfi.isValid(); ++mfi)
{
auto& ptile = ParticlesAt(lev, mfi);
auto& aos = ptile.GetArrayOfStructs();
auto pstruct_ptr = aos().dataPtr();
const size_t np = aos.numParticles();
const auto& ptile = ParticlesAt(lev, mfi);
const size_t np = ptile.numParticles();

const Box& box = mfi.validbox();

using index_type = typename decltype(m_bins)::index_type;
Gpu::DeviceVector<index_type> perm;
PermutationForDeposition<index_type>(perm, np, pstruct_ptr, box, geom, idx_type);
PermutationForDeposition<index_type>(perm, np, ptile, box, geom, idx_type);
ReorderParticles(lev, mfi, perm.dataPtr());
}
}
Expand Down
2 changes: 1 addition & 1 deletion Src/Particle/AMReX_ParticleTile.H
Original file line number Diff line number Diff line change
Expand Up @@ -316,7 +316,7 @@ struct ConstSoAParticle : SoAParticleBase
//functions to get positions of the particle in the SOA data

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
RealVect pos () const & {return RealVect(AMREX_D_DECL(this->m_constparticle_tile_data->m_rdata[0][m_index], this->m_constparticle_tile_data.m_rdata[1][m_index], this->m_constparticle_tile_data->m_rdata[2][m_index]));}
RealVect pos () const & {return RealVect(AMREX_D_DECL(this->m_constparticle_tile_data.m_rdata[0][m_index], this->m_constparticle_tile_data.m_rdata[1][m_index], this->m_constparticle_tile_data.m_rdata[2][m_index]));}

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
const RealType& pos (int position_index) const &
Expand Down
22 changes: 13 additions & 9 deletions Src/Particle/AMReX_ParticleUtil.H
Original file line number Diff line number Diff line change
Expand Up @@ -881,9 +881,9 @@ void PermutationForDeposition (Gpu::DeviceVector<index_type>& perm, index_type n
Gpu::Device::streamSynchronize();
}

template <class index_type, class pos_struct>
template <class index_type, class PTile>
void PermutationForDeposition (Gpu::DeviceVector<index_type>& perm, index_type nitems,
pos_struct const* v, Box bx, Geometry geom, const IntVect idx_type)
const PTile& ptile, Box bx, Geometry geom, const IntVect idx_type)
{
AMREX_ALWAYS_ASSERT(idx_type.allGE(IntVect(0)) && idx_type.allLE(IntVect(2)));

Expand All @@ -904,17 +904,21 @@ void PermutationForDeposition (Gpu::DeviceVector<index_type>& perm, index_type n
const int ref_product = AMREX_D_TERM(refine_vect[0], * refine_vect[1], * refine_vect[2]);
const IntVect ref_offset(AMREX_D_DECL(1, refine_vect[0], refine_vect[0] * refine_vect[1]));

auto ptd = ptile.getConstParticleTileData();
using ParticleType = typename PTile::ParticleType::ConstType;
PermutationForDeposition<index_type>(perm, nitems, bx.numPts() * ref_product,
[=] AMREX_GPU_DEVICE (index_type idx) noexcept
{
IntVect iv = ((v[idx].pos() - pos_offset) * dxi).round();
{
const auto& p = make_particle<ParticleType>{}(ptd,idx);

IntVect iv = ((p.pos() - pos_offset) * dxi).round();

IntVect iv_coarse = iv / refine_vect;
IntVect iv_remainder = iv - iv_coarse * refine_vect;
IntVect iv_coarse = iv / refine_vect;
IntVect iv_remainder = iv - iv_coarse * refine_vect;

iv_coarse = iv_coarse.max(bx.smallEnd());
iv_coarse = iv_coarse.min(bx.bigEnd());
return bx.index(iv_coarse) + bx.numPts() * (iv_remainder * ref_offset).sum();
iv_coarse = iv_coarse.max(bx.smallEnd());
iv_coarse = iv_coarse.min(bx.bigEnd());
return bx.index(iv_coarse) + bx.numPts() * (iv_remainder * ref_offset).sum();
});
}

Expand Down

0 comments on commit e7e2e1e

Please sign in to comment.