diff options
author | 2021-01-27 10:02:18 -0800 | |
---|---|---|
committer | 2021-01-27 10:02:18 -0800 | |
commit | 85bc4610adbdd5a48f1cbe666f11db6b72a781c0 (patch) | |
tree | 1616b2a7082196182a853ac85a336974f80dc680 /Source/Particles/WarpXParticleContainer.cpp | |
parent | 8c1c80cedbb86df828f49a3616681dc7d73c3946 (diff) | |
download | WarpX-85bc4610adbdd5a48f1cbe666f11db6b72a781c0.tar.gz WarpX-85bc4610adbdd5a48f1cbe666f11db6b72a781c0.tar.zst WarpX-85bc4610adbdd5a48f1cbe666f11db6b72a781c0.zip |
Take time step into account to compute guard cells for J and rho (#1607)
* Use IntVect for ng_J and ng_rho
* Compute guard cells for J and rho based on dt
* Reset some CI benchmarks
* Fix rebase commit
* Add back +1 cell for rho: fix remaining out-of-bound accesses
* Simplify ASSERTS using new interface of amrex::numParticlesOutOfRange
Diffstat (limited to 'Source/Particles/WarpXParticleContainer.cpp')
-rw-r--r-- | Source/Particles/WarpXParticleContainer.cpp | 87 |
1 files changed, 50 insertions, 37 deletions
diff --git a/Source/Particles/WarpXParticleContainer.cpp b/Source/Particles/WarpXParticleContainer.cpp index b0bbb75ef..1d17e534e 100644 --- a/Source/Particles/WarpXParticleContainer.cpp +++ b/Source/Particles/WarpXParticleContainer.cpp @@ -249,30 +249,37 @@ WarpXParticleContainer::DepositCurrent(WarpXParIter& pti, // Number of guard cells for local deposition of J WarpX& warpx = WarpX::GetInstance(); - const int ng_J = warpx.get_ng_depos_J().max(); - - // Extract deposition order (same order along all directions) and check that - // particles shape fits within the guard cells. - // NOTE: In specific situations where the staggering of J and the current - // deposition algorithm are not trivial, this check might be too relaxed - // and we might include a particle that should deposit part of its current - // in a neighboring box. However, this should catch particles traveling many - // cells away, for example with algorithms that allow for large time steps. - const int shape_extent = static_cast<int>(WarpX::nox / 2); + + const amrex::IntVect& ng_J = warpx.get_ng_depos_J(); + + // Extract deposition order and check that particles shape fits within the guard cells. + // NOTE: In specific situations where the staggering of J and the current deposition algorithm + // are not trivial, this check might be too relaxed and we might include a particle that should + // deposit part of its current in a neighboring box. However, this should catch particles + // traveling many cells away, for example with algorithms that allow for large time steps. + +#if (AMREX_SPACEDIM == 2) + const amrex::IntVect shape_extent = amrex::IntVect(static_cast<int>(WarpX::nox/2), + static_cast<int>(WarpX::noz/2)); +#elif (AMREX_SPACEDIM == 3) + const amrex::IntVect shape_extent = amrex::IntVect(static_cast<int>(WarpX::nox/2), + static_cast<int>(WarpX::noy/2), + static_cast<int>(WarpX::noz/2)); +#endif + + // On CPU: particles deposit on tile arrays, which have a small number of guard cells ng_J + // On GPU: particles deposit directly on the J arrays, which usually have a larger number of guard cells #ifndef AMREX_USE_GPU - // On CPU: Particles deposit on tile arrays, - // which have a small number of guard cells ng_J - AMREX_ALWAYS_ASSERT_WITH_MESSAGE( - amrex::numParticlesOutOfRange(pti, ng_J - shape_extent) == 0, - "Particles shape does not fit within tile used for local current deposition"); + const amrex::IntVect range = ng_J - shape_extent; #else - // On GPU: Particles deposit directly on the J arrays, - // which usually have a larger number of guard cells - AMREX_ALWAYS_ASSERT_WITH_MESSAGE( - amrex::numParticlesOutOfRange(pti, jx->nGrowVect().min() - shape_extent) == 0, - "Particles shape does not fit within guard cells used for current deposition"); + // Jx, Jy and Jz have the same number of guard cells, hence it is sufficient to check for Jx + const amrex::IntVect range = jx->nGrowVect() - shape_extent; #endif + AMREX_ALWAYS_ASSERT_WITH_MESSAGE( + amrex::numParticlesOutOfRange(pti, range) == 0, + "Particles shape does not fit within tile (CPU) or guard cells (GPU) used for current deposition"); + const std::array<Real,3>& dx = WarpX::CellSize(std::max(depos_lev,0)); Real q = this->charge; @@ -466,28 +473,34 @@ WarpXParticleContainer::DepositCharge (WarpXParIter& pti, RealVector& wp, // Number of guard cells for local deposition of rho WarpX& warpx = WarpX::GetInstance(); - const int ng_rho = warpx.get_ng_depos_rho().max(); + const amrex::IntVect& ng_rho = warpx.get_ng_depos_rho(); + + // Extract deposition order and check that particles shape fits within the guard cells. + // NOTE: In specific situations where the staggering of rho and the charge deposition algorithm + // are not trivial, this check might be too strict and we might need to relax it, as currently + // done for the current deposition. + +#if (AMREX_SPACEDIM == 2) + const amrex::IntVect shape_extent = amrex::IntVect(static_cast<int>(WarpX::nox/2+1), + static_cast<int>(WarpX::noz/2+1)); +#elif (AMREX_SPACEDIM == 3) + const amrex::IntVect shape_extent = amrex::IntVect(static_cast<int>(WarpX::nox/2+1), + static_cast<int>(WarpX::noy/2+1), + static_cast<int>(WarpX::noz/2+1)); +#endif - // Extract deposition order (same order along all directions) and check that - // particles shape fits within the guard cells. - // NOTE: In specific situations where the staggering of rho and the charge - // deposition algorithm are not trivial, this check might be too strict and - // we might need to relax it, as currently done for the current deposition. - const int shape_extent = static_cast<int>(WarpX::nox / 2 + 1); + // On CPU: particles deposit on tile arrays, which have a small number of guard cells ng_rho + // On GPU: particles deposit directly on the rho array, which usually have a larger number of guard cells #ifndef AMREX_USE_GPU - // On CPU: Particles deposit on tile arrays, - // which have a small number of guard cells - AMREX_ALWAYS_ASSERT_WITH_MESSAGE( - amrex::numParticlesOutOfRange(pti, ng_rho - shape_extent) == 0, - "Particles shape does not fit within tile used for local charge deposition"); + const amrex::IntVect range = ng_rho - shape_extent; #else - // On GPU: Particles deposit directly on the rho arrays, - // which usually have a larger number of guard cells - AMREX_ALWAYS_ASSERT_WITH_MESSAGE( - amrex::numParticlesOutOfRange(pti, rho->nGrowVect().min() - shape_extent) == 0, - "Particles shape does not fit within guard cells used for charge deposition"); + const amrex::IntVect range = rho->nGrowVect() - shape_extent; #endif + AMREX_ALWAYS_ASSERT_WITH_MESSAGE( + amrex::numParticlesOutOfRange(pti, range) == 0, + "Particles shape does not fit within tile (CPU) or guard cells (GPU) used for charge deposition"); + const std::array<Real,3>& dx = WarpX::CellSize(std::max(depos_lev,0)); const Real q = this->charge; |