WarpX
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>
64  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
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 
101  QuantumSynchrotronGetOpticalDepth opt_depth_functor,
102  int const opt_depth_runtime_comp,
103  QuantumSynchrotronPhotonEmission const emission_functor,
104  const WarpXParIter& a_pti, int lev, amrex::IntVect ngE,
105  amrex::FArrayBox const& exfab,
106  amrex::FArrayBox const& eyfab,
107  amrex::FArrayBox const& ezfab,
108  amrex::FArrayBox const& bxfab,
109  amrex::FArrayBox const& byfab,
110  amrex::FArrayBox const& bzfab,
111  amrex::Array<amrex::Real,3> v_galilean,
112  int a_offset = 0);
113 
124  template <typename DstData, typename SrcData>
125  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
126  void operator() (DstData& dst, SrcData& src, int i_src, int i_dst,
127  amrex::RandomEngine const& engine) const noexcept
128  {
129  using namespace amrex;
130 
131  // gather E and B
132  amrex::ParticleReal xp, yp, zp;
133  m_get_position(i_src, xp, yp, zp);
134 
135  amrex::ParticleReal ex = 0._rt, ey = 0._rt, ez = 0._rt;
136  m_get_externalE(i_src, ex, ey, ez);
137 
138  amrex::ParticleReal bx = 0._rt, by = 0._rt, bz = 0._rt;
139  m_get_externalB(i_src, bx, by, bz);
140 
141  doGatherShapeN(xp, yp, zp, ex, ey, ez, bx, by, bz,
142  m_ex_arr, m_ey_arr, m_ez_arr, m_bx_arr, m_by_arr, m_bz_arr,
143  m_ex_type, m_ey_type, m_ez_type, m_bx_type, m_by_type, m_bz_type,
144  m_dx_arr, m_xyzmin_arr, m_lo, m_n_rz_azimuthal_modes,
145  m_nox, m_galerkin_interpolation);
146 
147  auto& ux = src.m_rdata[PIdx::ux][i_src];
148  auto& uy = src.m_rdata[PIdx::uy][i_src];
149  auto& uz = src.m_rdata[PIdx::uz][i_src];
150  auto& g_ux = dst.m_rdata[PIdx::ux][i_dst];
151  auto& g_uy = dst.m_rdata[PIdx::uy][i_dst];
152  auto& g_uz = dst.m_rdata[PIdx::uz][i_dst];
153  m_emission_functor(
154  ux, uy, uz,
155  ex, ey, ez,
156  bx, by, bz,
157  g_ux, g_uy, g_uz,
158  engine);
159 
160  //Initialize the optical depth component of the source species.
161  src.m_runtime_rdata[m_opt_depth_runtime_comp][i_src] =
162  m_opt_depth_functor(engine);
163  }
164 
165 private:
169  const int m_opt_depth_runtime_comp = 0;
177 
178  amrex::Array4<const amrex::Real> m_ex_arr;
179  amrex::Array4<const amrex::Real> m_ey_arr;
180  amrex::Array4<const amrex::Real> m_ez_arr;
181  amrex::Array4<const amrex::Real> m_bx_arr;
182  amrex::Array4<const amrex::Real> m_by_arr;
183  amrex::Array4<const amrex::Real> m_bz_arr;
184 
185  amrex::IndexType m_ex_type;
186  amrex::IndexType m_ey_type;
187  amrex::IndexType m_ez_type;
188  amrex::IndexType m_bx_type;
189  amrex::IndexType m_by_type;
190  amrex::IndexType m_bz_type;
191 
192  amrex::GpuArray<amrex::Real, 3> m_dx_arr;
193  amrex::GpuArray<amrex::Real, 3> m_xyzmin_arr;
194 
196  int m_nox;
198 
199  amrex::Dim3 m_lo;
200 };
201 
202 
215 template <typename PTile>
217  PTile& ptile,
218  const int old_size, const int num_added,
219  const amrex::ParticleReal energy_threshold)
220 {
221  auto pp = ptile.GetArrayOfStructs()().data() + old_size;
222 
223  const auto& soa = ptile.GetStructOfArrays();
224  const auto p_ux = soa.GetRealData(PIdx::ux).data() + old_size;
225  const auto p_uy = soa.GetRealData(PIdx::uy).data() + old_size;
226  const auto p_uz = soa.GetRealData(PIdx::uz).data() + old_size;
227 
228  //The square of the energy threshold
229  const auto energy_threshold2 = std::max(
230  energy_threshold*energy_threshold,
231  std::numeric_limits<amrex::ParticleReal>::min());
232 
233  amrex::ParallelFor(num_added, [=] AMREX_GPU_DEVICE (int ip) noexcept
234  {
235  auto& p = pp[ip];
236 
237  const auto ux = p_ux[ip];
238  const auto uy = p_uy[ip];
239  const auto uz = p_uz[ip];
240 
241  //The square of the photon energy (in SI units)
242  // ( Particle momentum is stored as gamma * velocity.)
243  constexpr amrex::ParticleReal me_c = PhysConst::m_e*PhysConst::c;
244  const auto phot_energy2 = (ux*ux + uy*uy + uz*uz)*me_c*me_c;
245 
246  if (phot_energy2 < energy_threshold2){
247  p.id() = - 1;
248  }
249  });
250 }
251 
252 
253 #endif //QED_PHOTON_EMISSION_H_
GetExternalEField m_get_externalE
Definition: QEDPhotonEmission.H:175
int m_nox
Definition: QEDPhotonEmission.H:196
Definition: WarpXParticleContainer_fwd.H:28
int m_opt_depth_runtime_comp
Definition: QEDPhotonEmission.H:75
GetExternalBField m_get_externalB
Definition: QEDPhotonEmission.H:176
amrex::Array4< const amrex::Real > m_ez_arr
Definition: QEDPhotonEmission.H:180
data
Definition: run_alltests_1node.py:320
amrex::IndexType m_ex_type
Definition: QEDPhotonEmission.H:185
Definition: QuantumSyncEngineWrapper.H:77
amrex::GpuArray< amrex::Real, 3 > m_xyzmin_arr
Definition: QEDPhotonEmission.H:193
Definition: WarpXParticleContainer_fwd.H:28
bool m_galerkin_interpolation
Definition: QEDPhotonEmission.H:195
def uz
Definition: read_lab_particles.py:29
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
amrex::IndexType m_bz_type
Definition: QEDPhotonEmission.H:190
GetParticlePosition m_get_position
Definition: QEDPhotonEmission.H:174
amrex::Array4< const amrex::Real > m_ex_arr
Definition: QEDPhotonEmission.H:178
const QuantumSynchrotronPhotonEmission m_emission_functor
Definition: QEDPhotonEmission.H:172
Functor that can be used to assign the external B field to a particle inside a ParallelFor kernel...
Definition: GetExternalFields.H:159
Definition: WarpXParticleContainer_fwd.H:28
amrex::Array4< const amrex::Real > m_bx_arr
Definition: QEDPhotonEmission.H:181
amrex::Array4< const amrex::Real > m_bz_arr
Definition: QEDPhotonEmission.H:183
i
Definition: check_interp_points_and_weights.py:171
amrex::IndexType m_ey_type
Definition: QEDPhotonEmission.H:186
amrex::GpuArray< amrex::Real, 3 > m_dx_arr
Definition: QEDPhotonEmission.H:192
amrex::Array4< const amrex::Real > m_ey_arr
Definition: QEDPhotonEmission.H:179
amrex::Array4< const amrex::Real > m_by_arr
Definition: QEDPhotonEmission.H:182
int m_n_rz_azimuthal_modes
Definition: QEDPhotonEmission.H:197
Transform functor for the QED photon emission process.
Definition: QEDPhotonEmission.H:81
amrex::IndexType m_bx_type
Definition: QEDPhotonEmission.H:188
amrex::IndexType m_ez_type
Definition: QEDPhotonEmission.H:187
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
Definition: QuantumSyncEngineWrapper.H:180
amrex::IndexType m_by_type
Definition: QEDPhotonEmission.H:189
PhotonEmissionFilterFunc(int const opt_depth_runtime_comp)
Constructor of the PhotonEmissionFilterFunc functor.
Definition: QEDPhotonEmission.H:51
Functor that can be used to extract the positions of the macroparticles inside a ParallelFor kernel...
Definition: GetAndSetPosition.H:48
Filter functor for the QED photon emission process.
Definition: QEDPhotonEmission.H:42
const QuantumSynchrotronGetOpticalDepth m_opt_depth_functor
Definition: QEDPhotonEmission.H:167
def ux
Definition: read_lab_particles.py:28
amrex::Dim3 m_lo
Definition: QEDPhotonEmission.H:199
Definition: WarpXParticleContainer.H:58
Functor that can be used to assign the external E field to a particle inside a ParallelFor kernel...
Definition: GetExternalFields.H:149
Definition: BreitWheelerEngineWrapper.H:35
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:216