WarpX
|
#include "EmbeddedBoundary/DistanceToEB.H"
#include "Particles/Pusher/GetAndSetPosition.H"
#include <ablastr/particles/NodalFieldGather.H>
#include <AMReX.H>
#include <AMReX_MultiFab.H>
#include <AMReX_Particle.H>
#include <AMReX_RandomEngine.H>
#include <AMReX_REAL.H>
#include <AMReX_TypeTraits.H>
#include <AMReX_Vector.H>
#include <type_traits>
#include <utility>
Go to the source code of this file.
Functions | |
template<class PC , class F , std::enable_if_t< amrex::IsParticleContainer< PC >::value, int > foo = 0> | |
void | scrapeParticlesAtEB (PC &pc, const amrex::Vector< const amrex::MultiFab * > &distance_to_eb, int lev, F &&f) |
Interact particles with the embedded boundary walls. More... | |
template<class PC , class F , std::enable_if_t< amrex::IsParticleContainer< PC >::value, int > foo = 0> | |
void | scrapeParticlesAtEB (PC &pc, const amrex::Vector< const amrex::MultiFab * > &distance_to_eb, F &&f) |
Interact particles with the embedded boundary walls. More... | |
template<class PC , class F , std::enable_if_t< amrex::IsParticleContainer< PC >::value, int > foo = 0> | |
void | scrapeParticlesAtEB (PC &pc, const amrex::Vector< const amrex::MultiFab * > &distance_to_eb, int lev_min, int lev_max, F &&f) |
Interact particles with the embedded boundary walls. More... | |
void scrapeParticlesAtEB | ( | PC & | pc, |
const amrex::Vector< const amrex::MultiFab * > & | distance_to_eb, | ||
F && | f | ||
) |
Interact particles with the embedded boundary walls.
This function detects which particles have entered into the region covered by the embedded boundaries and applies an operation on those that have. Boundary collision detection is performed using a signed distance function, which is generated automatically when WarpX is compiled with EB support.
The operation to be performed is specified by the callable function passed in to this function as an argument. This function can access the position at which the particle hit the boundary, and also the associated normal vector. Particles can be absorbed
by setting their ids to negative to flag them for removal. Likewise, they can be reflected back into the domain by modifying their data appropriately and leaving their ids alone.
This version operates over all the levels in the pc.
pc | a type of amrex ParticleContainer |
F | a callable type, e.g. a lambda function or functor |
pc | the particle container to test for boundary interactions. |
distance_to_eb | a set of MultiFabs that store the signed distance function |
f | the callable that defines what to do when a particle hits the boundary. The form of the callable should model: template <typename PData> AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void operator() (const PData& ptd, int i, const amrex::RealVect& pos, const amrex::RealVect& normal, amrex::RandomEngine const& engine); where ptd is the particle tile, i the index of the particle operated on, pos and normal the location of the collision and the boundary normal vector. engine is for random number generation, if needed. |
void scrapeParticlesAtEB | ( | PC & | pc, |
const amrex::Vector< const amrex::MultiFab * > & | distance_to_eb, | ||
int | lev, | ||
F && | f | ||
) |
Interact particles with the embedded boundary walls.
This function detects which particles have entered into the region covered by the embedded boundaries and applies an operation on those that have. Boundary collision detection is performed using a signed distance function, which is generated automatically when WarpX is compiled with EB support.
The operation to be performed is specified by the callable function passed in to this function as an argument. This function can access the position at which the particle hit the boundary, and also the associated normal vector. Particles can be absorbed
by setting their ids to negative to flag them for removal. Likewise, they can be reflected back into the domain by modifying their data appropriately and leaving their ids alone.
This version operates only at the specified level.
pc | a type of amrex ParticleContainer |
F | a callable type, e.g. a lambda function or functor |
pc | the particle container to test for boundary interactions. |
distance_to_eb | a set of MultiFabs that store the signed distance function |
lev | the mesh refinement level to work on. |
f | the callable that defines what to do when a particle hits the boundary. The form of the callable should model: template <typename PData> AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void operator() (const PData& ptd, int i, const amrex::RealVect& pos, const amrex::RealVect& normal, amrex::RandomEngine const& engine); where ptd is the particle tile, i the index of the particle operated on, pos and normal the location of the collision and the boundary normal vector. engine is for random number generation, if needed. |
void scrapeParticlesAtEB | ( | PC & | pc, |
const amrex::Vector< const amrex::MultiFab * > & | distance_to_eb, | ||
int | lev_min, | ||
int | lev_max, | ||
F && | f | ||
) |
Interact particles with the embedded boundary walls.
This function detects which particles have entered into the region covered by the embedded boundaries and applies an operation on those that have. Boundary collision detection is performed using a signed distance function, which is generated automatically when WarpX is compiled with EB support.
The operation to be performed is specified by the callable function passed in to this function as an argument. This function can access the position at which the particle hit the boundary, and also the associated normal vector. Particles can be absorbed
by setting their ids to negative to flag them for removal. Likewise, the can be reflected back into the domain by modifying their data appropriately and leaving their ids alone.
This version operates only at the specified levels.
pc | a type of amrex ParticleContainer |
F | a callable type, e.g. a lambda function or functor |
pc | the particle container to test for boundary interactions. |
distance_to_eb | a set of MultiFabs that store the signed distance function |
lev_min | the minimum mesh refinement level to work on. |
lev_max | the maximum mesh refinement level to work on. |
f | the callable that defines what to do when a particle hits the boundary. The form of the callable should model: template <typename PData> AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void operator() (const PData& ptd, int i, const amrex::RealVect& pos, const amrex::RealVect& normal, amrex::RandomEngine const& engine); where ptd is the particle tile, i the index of the particle operated on, pos and normal the location of the collision and the boundary normal vector. engine is for random number generation, if needed. |