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>
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 opt_depth_runtime_comp,
113  QuantumSynchrotronPhotonEmission 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  amrex::Vector<amrex::ParticleReal>& E_external_particle,
122  amrex::Vector<amrex::ParticleReal>& B_external_particle,
123  int a_offset = 0);
124 
135  template <typename DstData, typename SrcData>
137  void operator() (DstData& dst, SrcData& src, int i_src, int i_dst,
138  amrex::RandomEngine const& engine) const noexcept
139  {
140  using namespace amrex;
141 
142  // gather E and B
143  amrex::ParticleReal xp, yp, zp;
144  m_get_position(i_src, xp, yp, zp);
145 
146  amrex::ParticleReal ex = m_Ex_external_particle;
147  amrex::ParticleReal ey = m_Ey_external_particle;
148  amrex::ParticleReal ez = m_Ez_external_particle;
149  amrex::ParticleReal bx = m_Bx_external_particle;
150  amrex::ParticleReal by = m_By_external_particle;
151  amrex::ParticleReal bz = m_Bz_external_particle;
152 
153  m_get_externalEB(i_src, ex, ey, ez, bx, by, bz);
154 
155  doGatherShapeN(xp, yp, zp, ex, ey, ez, bx, by, bz,
156  m_ex_arr, m_ey_arr, m_ez_arr, m_bx_arr, m_by_arr, m_bz_arr,
157  m_ex_type, m_ey_type, m_ez_type, m_bx_type, m_by_type, m_bz_type,
158  m_dx_arr, m_xyzmin_arr, m_lo, m_n_rz_azimuthal_modes,
159  m_nox, m_galerkin_interpolation);
160 
161  auto& ux = src.m_rdata[PIdx::ux][i_src];
162  auto& uy = src.m_rdata[PIdx::uy][i_src];
163  auto& uz = src.m_rdata[PIdx::uz][i_src];
164  auto& g_ux = dst.m_rdata[PIdx::ux][i_dst];
165  auto& g_uy = dst.m_rdata[PIdx::uy][i_dst];
166  auto& g_uz = dst.m_rdata[PIdx::uz][i_dst];
167  m_emission_functor(
168  ux, uy, uz,
169  ex, ey, ez,
170  bx, by, bz,
171  g_ux, g_uy, g_uz,
172  engine);
173 
174  //Initialize the optical depth component of the source species.
175  src.m_runtime_rdata[m_opt_depth_runtime_comp][i_src] =
176  m_opt_depth_functor(engine);
177  }
178 
179 private:
183  const int m_opt_depth_runtime_comp = 0;
190  amrex::ParticleReal m_Ex_external_particle;
191  amrex::ParticleReal m_Ey_external_particle;
192  amrex::ParticleReal m_Ez_external_particle;
193  amrex::ParticleReal m_Bx_external_particle;
194  amrex::ParticleReal m_By_external_particle;
195  amrex::ParticleReal m_Bz_external_particle;
196 
203 
210 
213 
215  int m_nox;
217 
219 };
220 
221 
234 template <typename PTile>
236  PTile& ptile,
237  const int old_size, const int num_added,
238  const amrex::ParticleReal energy_threshold)
239 {
240  auto pp = ptile.GetArrayOfStructs()().data() + old_size;
241 
242  const auto& soa = ptile.GetStructOfArrays();
243  const auto p_ux = soa.GetRealData(PIdx::ux).data() + old_size;
244  const auto p_uy = soa.GetRealData(PIdx::uy).data() + old_size;
245  const auto p_uz = soa.GetRealData(PIdx::uz).data() + old_size;
246 
247  //The square of the energy threshold
248  const auto energy_threshold2 = std::max(
249  energy_threshold*energy_threshold,
250  std::numeric_limits<amrex::ParticleReal>::min());
251 
252  amrex::ParallelFor(num_added, [=] AMREX_GPU_DEVICE (int ip) noexcept
253  {
254  auto& p = pp[ip];
255 
256  const auto ux = p_ux[ip];
257  const auto uy = p_uy[ip];
258  const auto uz = p_uz[ip];
259 
260  //The square of the photon energy (in SI units)
261  // ( Particle momentum is stored as gamma * velocity.)
262  constexpr amrex::ParticleReal me_c = PhysConst::m_e*PhysConst::c;
263  const auto phot_energy2 = (ux*ux + uy*uy + uz*uz)*me_c*me_c;
264 
265  if (phot_energy2 < energy_threshold2){
266  p.id() = - 1;
267  }
268  });
269 }
270 
271 
272 #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:235
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
GetParticlePosition< PIdx > m_get_position
Definition: QEDPhotonEmission.H:188
amrex::ParticleReal m_Ey_external_particle
Definition: QEDPhotonEmission.H:191
amrex::Array4< const amrex::Real > m_ey_arr
Definition: QEDPhotonEmission.H:198
amrex::GpuArray< amrex::Real, 3 > m_xyzmin_arr
Definition: QEDPhotonEmission.H:212
amrex::ParticleReal m_By_external_particle
Definition: QEDPhotonEmission.H:194
bool m_galerkin_interpolation
Definition: QEDPhotonEmission.H:214
amrex::IndexType m_bz_type
Definition: QEDPhotonEmission.H:209
int m_nox
Definition: QEDPhotonEmission.H:215
amrex::ParticleReal m_Ez_external_particle
Definition: QEDPhotonEmission.H:192
const QuantumSynchrotronGetOpticalDepth m_opt_depth_functor
Definition: QEDPhotonEmission.H:181
amrex::IndexType m_by_type
Definition: QEDPhotonEmission.H:208
amrex::Array4< const amrex::Real > m_ex_arr
Definition: QEDPhotonEmission.H:197
amrex::IndexType m_ex_type
Definition: QEDPhotonEmission.H:204
amrex::IndexType m_bx_type
Definition: QEDPhotonEmission.H:207
amrex::IndexType m_ey_type
Definition: QEDPhotonEmission.H:205
amrex::IndexType m_ez_type
Definition: QEDPhotonEmission.H:206
amrex::Array4< const amrex::Real > m_bx_arr
Definition: QEDPhotonEmission.H:200
amrex::GpuArray< amrex::Real, 3 > m_dx_arr
Definition: QEDPhotonEmission.H:211
int m_n_rz_azimuthal_modes
Definition: QEDPhotonEmission.H:216
amrex::Array4< const amrex::Real > m_ez_arr
Definition: QEDPhotonEmission.H:199
amrex::ParticleReal m_Bz_external_particle
Definition: QEDPhotonEmission.H:195
amrex::ParticleReal m_Ex_external_particle
Definition: QEDPhotonEmission.H:190
amrex::Array4< const amrex::Real > m_by_arr
Definition: QEDPhotonEmission.H:201
const QuantumSynchrotronPhotonEmission m_emission_functor
Definition: QEDPhotonEmission.H:186
GetExternalEBField m_get_externalEB
Definition: QEDPhotonEmission.H:189
amrex::ParticleReal m_Bx_external_particle
Definition: QEDPhotonEmission.H:193
amrex::Dim3 m_lo
Definition: QEDPhotonEmission.H:218
amrex::Array4< const amrex::Real > m_bz_arr
Definition: QEDPhotonEmission.H:202
Definition: QuantumSyncEngineWrapper.H:76
Definition: QuantumSyncEngineWrapper.H:179
Definition: WarpXParticleContainer.H:53
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< T >::value > ParallelFor(TypeList< CTOs... >, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
i
Definition: check_interp_points_and_weights.py:174
data
Definition: run_alltests_1node.py:325
Functor class that assigns external field values (E and B) to particles.
Definition: GetExternalFields.H:25
@ uz
Definition: NamedComponentParticleContainer.H:27
@ uy
Definition: NamedComponentParticleContainer.H:27
@ ux
Definition: NamedComponentParticleContainer.H:27