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 
11 #include "Utils/WarpXConst.H"
17 
18 #include <limits>
19 
31 {
32 public:
33 
39  PhotonEmissionFilterFunc(int const opt_depth_runtime_comp)
40  : m_opt_depth_runtime_comp(opt_depth_runtime_comp)
41  {}
42 
51  template <typename PData>
52  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
53  bool operator() (const PData& ptd, int const i, amrex::RandomEngine const&) const noexcept
54  {
55  using namespace amrex;
56 
57  const amrex::ParticleReal opt_depth =
58  ptd.m_runtime_rdata[m_opt_depth_runtime_comp][i];
59  return (opt_depth < 0.0_rt);
60  }
61 
62 private:
64 };
65 
70 {
71 
72 public:
73 
89  QuantumSynchrotronGetOpticalDepth opt_depth_functor,
90  int const opt_depth_runtime_comp,
91  QuantumSynchrotronPhotonEmission const emission_functor,
92  const WarpXParIter& a_pti, int lev, int ngE,
93  amrex::FArrayBox const& exfab,
94  amrex::FArrayBox const& eyfab,
95  amrex::FArrayBox const& ezfab,
96  amrex::FArrayBox const& bxfab,
97  amrex::FArrayBox const& byfab,
98  amrex::FArrayBox const& bzfab,
99  amrex::Array<amrex::Real,3> v_galilean,
100  int a_offset = 0);
101 
111  template <typename DstData, typename SrcData>
112  AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
113  void operator() (DstData& dst, SrcData& src, int i_src, int i_dst) const noexcept
114  {
115  using namespace amrex;
116 
117  // gather E and B
118  amrex::ParticleReal xp, yp, zp;
119  m_get_position(i_src, xp, yp, zp);
120 
121  amrex::ParticleReal ex = 0._rt, ey = 0._rt, ez = 0._rt;
122  m_get_externalE(i_src, ex, ey, ez);
123 
124  amrex::ParticleReal bx = 0._rt, by = 0._rt, bz = 0._rt;
125  m_get_externalB(i_src, bx, by, bz);
126 
127  doGatherShapeN(xp, yp, zp, ex, ey, ez, bx, by, bz,
128  m_ex_arr, m_ey_arr, m_ez_arr, m_bx_arr, m_by_arr, m_bz_arr,
129  m_ex_type, m_ey_type, m_ez_type, m_bx_type, m_by_type, m_bz_type,
130  m_dx_arr, m_xyzmin_arr, m_lo, m_n_rz_azimuthal_modes,
131  m_nox, m_galerkin_interpolation);
132 
133  auto& ux = src.m_rdata[PIdx::ux][i_src];
134  auto& uy = src.m_rdata[PIdx::uy][i_src];
135  auto& uz = src.m_rdata[PIdx::uz][i_src];
136  auto& g_ux = dst.m_rdata[PIdx::ux][i_dst];
137  auto& g_uy = dst.m_rdata[PIdx::uy][i_dst];
138  auto& g_uz = dst.m_rdata[PIdx::uz][i_dst];
139  m_emission_functor(
140  ux, uy, uz,
141  ex, ey, ez,
142  bx, by, bz,
143  g_ux, g_uy, g_uz);
144 
145  //Initialize the optical depth component of the source species.
146  src.m_runtime_rdata[m_opt_depth_runtime_comp][i_src] =
147  m_opt_depth_functor();
148  }
149 
150 private:
154  const int m_opt_depth_runtime_comp = 0;
162 
163  amrex::Array4<const amrex::Real> m_ex_arr;
164  amrex::Array4<const amrex::Real> m_ey_arr;
165  amrex::Array4<const amrex::Real> m_ez_arr;
166  amrex::Array4<const amrex::Real> m_bx_arr;
167  amrex::Array4<const amrex::Real> m_by_arr;
168  amrex::Array4<const amrex::Real> m_bz_arr;
169 
170  amrex::IndexType m_ex_type;
171  amrex::IndexType m_ey_type;
172  amrex::IndexType m_ez_type;
173  amrex::IndexType m_bx_type;
174  amrex::IndexType m_by_type;
175  amrex::IndexType m_bz_type;
176 
177  amrex::GpuArray<amrex::Real, 3> m_dx_arr;
178  amrex::GpuArray<amrex::Real, 3> m_xyzmin_arr;
179 
181  int m_nox;
183 
184  amrex::Dim3 m_lo;
185 };
186 
187 
200 template <typename PTile>
202  PTile& ptile,
203  const int old_size, const int num_added,
204  const amrex::ParticleReal energy_threshold)
205 {
206  auto pp = ptile.GetArrayOfStructs()().data() + old_size;
207 
208  const auto& soa = ptile.GetStructOfArrays();
209  const auto p_ux = soa.GetRealData(PIdx::ux).data() + old_size;
210  const auto p_uy = soa.GetRealData(PIdx::uy).data() + old_size;
211  const auto p_uz = soa.GetRealData(PIdx::uz).data() + old_size;
212 
213  //The square of the energy threshold
214  const auto energy_threshold2 = std::max(
215  energy_threshold*energy_threshold,
216  std::numeric_limits<amrex::ParticleReal>::min());
217 
218  amrex::ParallelFor(num_added, [=] AMREX_GPU_DEVICE (int ip) noexcept
219  {
220  auto& p = pp[ip];
221 
222  const auto ux = p_ux[ip];
223  const auto uy = p_uy[ip];
224  const auto uz = p_uz[ip];
225 
226  //The square of the photon energy (in SI units)
227  // ( Particle momentum is stored as gamma * velocity.)
228  constexpr amrex::ParticleReal me_c = PhysConst::m_e*PhysConst::c;
229  const auto phot_energy2 = (ux*ux + uy*uy + uz*uz)*me_c*me_c;
230 
231  if (phot_energy2 < energy_threshold2){
232  p.id() = - 1;
233  }
234  });
235 }
236 
237 
238 #endif //QED_PHOTON_EMISSION_H_
GetExternalEField m_get_externalE
Definition: QEDPhotonEmission.H:160
int m_nox
Definition: QEDPhotonEmission.H:181
int m_opt_depth_runtime_comp
Definition: QEDPhotonEmission.H:63
GetExternalBField m_get_externalB
Definition: QEDPhotonEmission.H:161
amrex::Array4< const amrex::Real > m_ez_arr
Definition: QEDPhotonEmission.H:165
data
Definition: run_alltests_1node.py:320
amrex::IndexType m_ex_type
Definition: QEDPhotonEmission.H:170
Definition: QuantumSyncEngineWrapper.H:65
amrex::GpuArray< amrex::Real, 3 > m_xyzmin_arr
Definition: QEDPhotonEmission.H:178
bool m_galerkin_interpolation
Definition: QEDPhotonEmission.H:180
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:53
amrex::IndexType m_bz_type
Definition: QEDPhotonEmission.H:175
GetParticlePosition m_get_position
Definition: QEDPhotonEmission.H:159
amrex::Array4< const amrex::Real > m_ex_arr
Definition: QEDPhotonEmission.H:163
const QuantumSynchrotronPhotonEmission m_emission_functor
Definition: QEDPhotonEmission.H:157
Functor that can be used to assign the external B field to a particle inside a ParallelFor kernel...
Definition: GetExternalFields.H:68
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 long n_rz_azimuthal_modes)
Field gather for a single particle.
Definition: FieldGather.H:38
amrex::Array4< const amrex::Real > m_bx_arr
Definition: QEDPhotonEmission.H:166
amrex::Array4< const amrex::Real > m_bz_arr
Definition: QEDPhotonEmission.H:168
i
Definition: check_interp_points_and_weights.py:171
amrex::IndexType m_ey_type
Definition: QEDPhotonEmission.H:171
Definition: WarpXParticleContainer.H:37
amrex::GpuArray< amrex::Real, 3 > m_dx_arr
Definition: QEDPhotonEmission.H:177
Definition: WarpXParticleContainer.H:37
amrex::Array4< const amrex::Real > m_ey_arr
Definition: QEDPhotonEmission.H:164
amrex::Array4< const amrex::Real > m_by_arr
Definition: QEDPhotonEmission.H:167
int m_n_rz_azimuthal_modes
Definition: QEDPhotonEmission.H:182
Transform functor for the QED photon emission process.
Definition: QEDPhotonEmission.H:69
amrex::IndexType m_bx_type
Definition: QEDPhotonEmission.H:173
Definition: WarpXParticleContainer.H:37
amrex::IndexType m_ez_type
Definition: QEDPhotonEmission.H:172
Definition: QuantumSyncEngineWrapper.H:168
amrex::IndexType m_by_type
Definition: QEDPhotonEmission.H:174
PhotonEmissionFilterFunc(int const opt_depth_runtime_comp)
Constructor of the PhotonEmissionFilterFunc functor.
Definition: QEDPhotonEmission.H:39
Functor that can be used to extract the positions of the macroparticles inside a ParallelFor kernel...
Definition: GetAndSetPosition.H:25
Filter functor for the QED photon emission process.
Definition: QEDPhotonEmission.H:30
const QuantumSynchrotronGetOpticalDepth m_opt_depth_functor
Definition: QEDPhotonEmission.H:152
def ux
Definition: read_lab_particles.py:28
amrex::Dim3 m_lo
Definition: QEDPhotonEmission.H:184
Definition: WarpXParticleContainer.H:76
Functor that can be used to assign the external E field to a particle inside a ParallelFor kernel...
Definition: GetExternalFields.H:58
Definition: PML.H:52
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:201