8#ifndef WARPX_QED_PHOTON_EMISSION_H_
9#define WARPX_QED_PHOTON_EMISSION_H_
64 template <
typename PData>
68 using namespace amrex;
70 const amrex::ParticleReal opt_depth =
72 return (opt_depth < 0.0_rt);
115 int opt_depth_runtime_comp,
138 template <
typename DstData,
typename SrcData>
140 void operator() (DstData& dst, SrcData& src,
int i_src,
int i_dst,
143 using namespace amrex;
146 amrex::ParticleReal xp, yp, zp;
164 auto& ux = src.m_rdata[
PIdx::ux][i_src];
165 auto& uy = src.m_rdata[
PIdx::uy][i_src];
166 auto& uz = src.m_rdata[
PIdx::uz][i_src];
167 auto& g_ux = dst.m_rdata[
PIdx::ux][i_dst];
168 auto& g_uy = dst.m_rdata[
PIdx::uy][i_dst];
169 auto& g_uz = dst.m_rdata[
PIdx::uz][i_dst];
237template <
typename PTile>
240 const int old_size,
const int num_added,
241 const amrex::ParticleReal energy_threshold)
243 auto& soa = ptile.GetStructOfArrays();
244 auto p_idcpu = soa.GetIdCPUData().data() + old_size;
245 const auto p_ux = soa.GetRealData(
PIdx::ux).data() + old_size;
246 const auto p_uy = soa.GetRealData(
PIdx::uy).data() + old_size;
247 const auto p_uz = soa.GetRealData(
PIdx::uz).data() + old_size;
250 const auto energy_threshold2 = std::max(
251 energy_threshold*energy_threshold,
252 std::numeric_limits<amrex::ParticleReal>::min());
256 const auto ux = p_ux[ip];
257 const auto uy = p_uy[ip];
258 const auto uz = p_uz[ip];
263 const auto phot_energy2 = (ux*ux + uy*uy + uz*uz)*me_c*me_c;
265 if (phot_energy2 < energy_threshold2) {
#define AMREX_FORCE_INLINE
#define AMREX_GPU_HOST_DEVICE
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void doGatherShapeN(const amrex::ParticleReal xp, const amrex::ParticleReal yp, const amrex::ParticleReal zp, amrex::ParticleReal &Exp, amrex::ParticleReal &Eyp, amrex::ParticleReal &Ezp, amrex::ParticleReal &Bxp, amrex::ParticleReal &Byp, amrex::ParticleReal &Bzp, amrex::Array4< amrex::Real const > const &ex_arr, amrex::Array4< amrex::Real const > const &ey_arr, amrex::Array4< amrex::Real const > const &ez_arr, amrex::Array4< amrex::Real const > const &bx_arr, amrex::Array4< amrex::Real const > const &by_arr, amrex::Array4< amrex::Real const > const &bz_arr, const amrex::IndexType ex_type, const amrex::IndexType ey_type, const amrex::IndexType ez_type, const amrex::IndexType bx_type, const amrex::IndexType by_type, const amrex::IndexType bz_type, const amrex::XDim3 &dinv, const amrex::XDim3 &xyzmin, const amrex::Dim3 &lo, const int n_rz_azimuthal_modes)
Field gather for a single particle.
Definition FieldGather.H:349
void cleanLowEnergyPhotons(PTile &ptile, const int old_size, const int num_added, const amrex::ParticleReal energy_threshold)
Free function to call to remove immediately low energy photons by setting their ID to -1....
Definition QEDPhotonEmission.H:238
int m_opt_depth_runtime_comp
Definition QEDPhotonEmission.H:76
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE bool operator()(const PData &ptd, int const i, amrex::RandomEngine const &) const noexcept
Functor call. This method determines if a given (electron or positron) particle should undergo QED ph...
Definition QEDPhotonEmission.H:66
PhotonEmissionFilterFunc(int const opt_depth_runtime_comp)
Constructor of the PhotonEmissionFilterFunc functor.
Definition QEDPhotonEmission.H:52
Definition QuantumSyncEngineWrapper.H:76
Definition QuantumSyncEngineWrapper.H:179
Definition WarpXParticleContainer.H:112
constexpr auto c
vacuum speed of light [m/s]
Definition constant.H:149
constexpr auto m_e
electron mass [kg]
Definition constant.H:161
constexpr std::uint64_t Invalid
std::enable_if_t< std::is_integral_v< T > > ParallelFor(TypeList< CTOs... > ctos, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
IndexTypeND< 3 > IndexType
Functor class that assigns external field values (E and B) to particles.
Definition GetExternalFields.H:23
Functor that can be used to extract the positions of the macroparticles inside a ParallelFor kernel.
Definition GetAndSetPosition.H:75
@ uz
Definition WarpXParticleContainer.H:70
@ uy
Definition WarpXParticleContainer.H:70
@ ux
Definition WarpXParticleContainer.H:70