ARTEMIS
QEDPhotonEmission.H
Go to the documentation of this file.
1 /* Copyright 2019 Luca Fedeli
2  *
3  * This file is part of WarpX.
4  *
5  * License: BSD-3-Clause-LBNL
6  */
7 
8 #ifndef QED_PHOTON_EMISSION_H_
9 #define QED_PHOTON_EMISSION_H_
10 
16 #include "Utils/WarpXConst.H"
17 
18 #include <AMReX_Array.H>
19 #include <AMReX_Array4.H>
20 #include <AMReX_Dim3.H>
21 #include <AMReX_Extension.H>
22 #include <AMReX_GpuLaunch.H>
23 #include <AMReX_GpuQualifiers.H>
24 #include <AMReX_IndexType.H>
25 #include <AMReX_REAL.H>
26 
27 #include <AMReX_BaseFwd.H>
28 
29 #include <algorithm>
30 #include <limits>
31 
43 {
44 public:
45 
51  PhotonEmissionFilterFunc(int const opt_depth_runtime_comp)
52  : m_opt_depth_runtime_comp(opt_depth_runtime_comp)
53  {}
54 
63  template <typename PData>
65  bool operator() (const PData& ptd, int const i, amrex::RandomEngine const&) const noexcept
66  {
67  using namespace amrex;
68 
69  const amrex::ParticleReal opt_depth =
70  ptd.m_runtime_rdata[m_opt_depth_runtime_comp][i];
71  return (opt_depth < 0.0_rt);
72  }
73 
74 private:
76 };
77 
82 {
83 
84 public:
85 
111  QuantumSynchrotronGetOpticalDepth opt_depth_functor,
112  int const opt_depth_runtime_comp,
113  QuantumSynchrotronPhotonEmission const emission_functor,
114  const WarpXParIter& a_pti, int lev, amrex::IntVect ngEB,
115  amrex::FArrayBox const& exfab,
116  amrex::FArrayBox const& eyfab,
117  amrex::FArrayBox const& ezfab,
118  amrex::FArrayBox const& bxfab,
119  amrex::FArrayBox const& byfab,
120  amrex::FArrayBox const& bzfab,
121  int a_offset = 0);
122 
133  template <typename DstData, typename SrcData>
135  void operator() (DstData& dst, SrcData& src, int i_src, int i_dst,
136  amrex::RandomEngine const& engine) const noexcept
137  {
138  using namespace amrex;
139 
140  // gather E and B
141  amrex::ParticleReal xp, yp, zp;
142  m_get_position(i_src, xp, yp, zp);
143 
144  amrex::ParticleReal ex = 0._rt, ey = 0._rt, ez = 0._rt;
145  amrex::ParticleReal bx = 0._rt, by = 0._rt, bz = 0._rt;
146  m_get_externalEB(i_src, ex, ey, ez, bx, by, bz);
147 
148  doGatherShapeN(xp, yp, zp, ex, ey, ez, bx, by, bz,
149  m_ex_arr, m_ey_arr, m_ez_arr, m_bx_arr, m_by_arr, m_bz_arr,
150  m_ex_type, m_ey_type, m_ez_type, m_bx_type, m_by_type, m_bz_type,
151  m_dx_arr, m_xyzmin_arr, m_lo, m_n_rz_azimuthal_modes,
152  m_nox, m_galerkin_interpolation);
153 
154  auto& ux = src.m_rdata[PIdx::ux][i_src];
155  auto& uy = src.m_rdata[PIdx::uy][i_src];
156  auto& uz = src.m_rdata[PIdx::uz][i_src];
157  auto& g_ux = dst.m_rdata[PIdx::ux][i_dst];
158  auto& g_uy = dst.m_rdata[PIdx::uy][i_dst];
159  auto& g_uz = dst.m_rdata[PIdx::uz][i_dst];
160  m_emission_functor(
161  ux, uy, uz,
162  ex, ey, ez,
163  bx, by, bz,
164  g_ux, g_uy, g_uz,
165  engine);
166 
167  //Initialize the optical depth component of the source species.
168  src.m_runtime_rdata[m_opt_depth_runtime_comp][i_src] =
169  m_opt_depth_functor(engine);
170  }
171 
172 private:
176  const int m_opt_depth_runtime_comp = 0;
183 
190 
197 
200 
202  int m_nox;
204 
206 };
207 
208 
221 template <typename PTile>
223  PTile& ptile,
224  const int old_size, const int num_added,
225  const amrex::ParticleReal energy_threshold)
226 {
227  auto pp = ptile.GetArrayOfStructs()().data() + old_size;
228 
229  const auto& soa = ptile.GetStructOfArrays();
230  const auto p_ux = soa.GetRealData(PIdx::ux).data() + old_size;
231  const auto p_uy = soa.GetRealData(PIdx::uy).data() + old_size;
232  const auto p_uz = soa.GetRealData(PIdx::uz).data() + old_size;
233 
234  //The square of the energy threshold
235  const auto energy_threshold2 = std::max(
236  energy_threshold*energy_threshold,
237  std::numeric_limits<amrex::ParticleReal>::min());
238 
239  amrex::ParallelFor(num_added, [=] AMREX_GPU_DEVICE (int ip) noexcept
240  {
241  auto& p = pp[ip];
242 
243  const auto ux = p_ux[ip];
244  const auto uy = p_uy[ip];
245  const auto uz = p_uz[ip];
246 
247  //The square of the photon energy (in SI units)
248  // ( Particle momentum is stored as gamma * velocity.)
249  constexpr amrex::ParticleReal me_c = PhysConst::m_e*PhysConst::c;
250  const auto phot_energy2 = (ux*ux + uy*uy + uz*uz)*me_c*me_c;
251 
252  if (phot_energy2 < energy_threshold2){
253  p.id() = - 1;
254  }
255  });
256 }
257 
258 
259 #endif //QED_PHOTON_EMISSION_H_
#define AMREX_FORCE_INLINE
#define AMREX_GPU_DEVICE
#define AMREX_GPU_HOST_DEVICE
amrex::ParmParse pp
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::GpuArray< amrex::Real, 3 > &dx, const amrex::GpuArray< amrex::Real, 3 > &xyzmin, const amrex::Dim3 &lo, const int n_rz_azimuthal_modes)
Field gather for a single particle.
Definition: FieldGather.H:38
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:222
Filter functor for the QED photon emission process.
Definition: QEDPhotonEmission.H:43
int m_opt_depth_runtime_comp
Definition: QEDPhotonEmission.H:75
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:65
PhotonEmissionFilterFunc(int const opt_depth_runtime_comp)
Constructor of the PhotonEmissionFilterFunc functor.
Definition: QEDPhotonEmission.H:51
Transform functor for the QED photon emission process.
Definition: QEDPhotonEmission.H:82
amrex::Array4< const amrex::Real > m_ey_arr
Definition: QEDPhotonEmission.H:185
amrex::GpuArray< amrex::Real, 3 > m_xyzmin_arr
Definition: QEDPhotonEmission.H:199
bool m_galerkin_interpolation
Definition: QEDPhotonEmission.H:201
amrex::IndexType m_bz_type
Definition: QEDPhotonEmission.H:196
int m_nox
Definition: QEDPhotonEmission.H:202
const QuantumSynchrotronGetOpticalDepth m_opt_depth_functor
Definition: QEDPhotonEmission.H:174
amrex::IndexType m_by_type
Definition: QEDPhotonEmission.H:195
amrex::Array4< const amrex::Real > m_ex_arr
Definition: QEDPhotonEmission.H:184
amrex::IndexType m_ex_type
Definition: QEDPhotonEmission.H:191
amrex::IndexType m_bx_type
Definition: QEDPhotonEmission.H:194
amrex::IndexType m_ey_type
Definition: QEDPhotonEmission.H:192
amrex::IndexType m_ez_type
Definition: QEDPhotonEmission.H:193
amrex::Array4< const amrex::Real > m_bx_arr
Definition: QEDPhotonEmission.H:187
amrex::GpuArray< amrex::Real, 3 > m_dx_arr
Definition: QEDPhotonEmission.H:198
int m_n_rz_azimuthal_modes
Definition: QEDPhotonEmission.H:203
amrex::Array4< const amrex::Real > m_ez_arr
Definition: QEDPhotonEmission.H:186
GetParticlePosition m_get_position
Definition: QEDPhotonEmission.H:181
amrex::Array4< const amrex::Real > m_by_arr
Definition: QEDPhotonEmission.H:188
const QuantumSynchrotronPhotonEmission m_emission_functor
Definition: QEDPhotonEmission.H:179
GetExternalEBField m_get_externalEB
Definition: QEDPhotonEmission.H:182
amrex::Dim3 m_lo
Definition: QEDPhotonEmission.H:205
amrex::Array4< const amrex::Real > m_bz_arr
Definition: QEDPhotonEmission.H:189
Definition: QuantumSyncEngineWrapper.H:76
Definition: QuantumSyncEngineWrapper.H:178
Definition: WarpXParticleContainer.H:52
list data
open text file and read data #####
Definition: Excitation_Flag_Generator.py:24
static constexpr auto c
vacuum speed of light [m/s]
Definition: constant.H:44
static constexpr auto m_e
electron mass [kg]
Definition: constant.H:52
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)
i
Definition: check_interp_points_and_weights.py:174
Functor class that assigns external field values (E and B) to particles.
Definition: GetExternalFields.H:25
Functor that can be used to extract the positions of the macroparticles inside a ParallelFor kernel.
Definition: GetAndSetPosition.H:53
@ uz
Definition: NamedComponentParticleContainer.H:25
@ uy
Definition: NamedComponentParticleContainer.H:25
@ ux
Definition: NamedComponentParticleContainer.H:25