ARTEMIS
MultiParticleContainer.H
Go to the documentation of this file.
1 /* Copyright 2019-2020 Andrew Myers, Ann Almgren, Axel Huebl
2  * David Grote, Jean-Luc Vay, Junmin Gu
3  * Luca Fedeli, Mathieu Lobet, Maxence Thevenet
4  * Remi Lehe, Revathi Jambunathan, Weiqun Zhang
5  * Yinjian Zhao
6  *
7  * This file is part of WarpX.
8  *
9  * License: BSD-3-Clause-LBNL
10  */
11 #ifndef WARPX_ParticleContainer_H_
12 #define WARPX_ParticleContainer_H_
13 
15 
16 #include "Evolve/WarpXDtType.H"
18 #ifdef WARPX_QED
21 #endif
23 #include "Utils/TextMsg.H"
24 #include "Utils/WarpXConst.H"
25 #include "WarpXParticleContainer.H"
26 #include "ParticleBoundaries.H"
27 
28 #include <AMReX_BLassert.H>
29 #include <AMReX_Box.H>
30 #include <AMReX_Config.H>
31 #include <AMReX_GpuControl.H>
32 #include <AMReX_INT.H>
33 #include <AMReX_MFIter.H>
34 #include <AMReX_REAL.H>
35 #include <AMReX_RealBox.H>
36 #include <AMReX_Vector.H>
37 
38 #include <AMReX_BaseFwd.H>
39 #include <AMReX_AmrCoreFwd.H>
40 
41 #include <algorithm>
42 #include <array>
43 #include <iosfwd>
44 #include <iterator>
45 #include <limits>
46 #include <memory>
47 #include <string>
48 #include <vector>
49 
65 {
66 
67 public:
68 
70 
72 
74  GetParticleContainer (int ispecies) const {return *allcontainers[ispecies];}
75 
77  GetParticleContainerPtr (int ispecies) const {return allcontainers[ispecies].get();}
78 
80  GetParticleContainerFromName (const std::string& name) const;
81 
82 #ifdef WARPX_USE_OPENPMD
83  std::unique_ptr<WarpXParticleContainer>& GetUniqueContainer(int ispecies) {
84  return allcontainers[ispecies];
85  }
86 #endif
87  std::array<amrex::ParticleReal, 3> meanParticleVelocity(int ispecies) {
88  return allcontainers[ispecies]->meanParticleVelocity();
89  }
90 
91  void AllocData ();
92 
93  void InitData ();
94 
96 
102  void Evolve (int lev,
103  const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
104  const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz,
107  amrex::MultiFab* rho, amrex::MultiFab* crho,
108  const amrex::MultiFab* cEx, const amrex::MultiFab* cEy, const amrex::MultiFab* cEz,
109  const amrex::MultiFab* cBx, const amrex::MultiFab* cBy, const amrex::MultiFab* cBz,
110  amrex::Real t, amrex::Real dt, DtType a_dt_type=DtType::Full, bool skip_deposition=false);
111 
117  void PushX (amrex::Real dt);
118 
125  void PushP (int lev, amrex::Real dt,
126  const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
127  const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz);
128 
135  std::unique_ptr<amrex::MultiFab> GetZeroChargeDensity(const int lev);
136 
146  void
147  DepositCharge (amrex::Vector<std::unique_ptr<amrex::MultiFab> >& rho,
148  const amrex::Real relative_time);
149 
161  void
162  DepositCurrent (amrex::Vector<std::array< std::unique_ptr<amrex::MultiFab>, 3 > >& J,
163  const amrex::Real dt, const amrex::Real relative_time);
164 
168 
169  std::unique_ptr<amrex::MultiFab> GetChargeDensity(int lev, bool local = false);
170 
171  void doFieldIonization (int lev,
172  const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
173  const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz);
174 
175  //void doCollisions (amrex::Real cur_time, amrex::Real dt);
176 
182  void doResampling (const int timestep);
183 
184 #ifdef WARPX_QED
192  void doQEDSchwinger ();
193 
199 #endif
200 
201  void Restart (const std::string& dir);
202 
203  void PostRestart ();
204 
205  void ReadHeader (std::istream& is);
206 
207  void WriteHeader (std::ostream& os) const;
208 
209  void SortParticlesByBin (amrex::IntVect bin_size);
210 
211  void Redistribute ();
212 
213  void defineAllParticleTiles ();
214 
215  void RedistributeLocal (const int num_ghost);
216 
219  void ApplyBoundaryConditions ();
220 
229 
231 
232  void Increment (amrex::MultiFab& mf, int lev);
233 
234  void SetParticleBoxArray (int lev, amrex::BoxArray& new_ba);
236 
237  int nSpecies() const {return species_names.size();}
238 
243  void SetDoBackTransformedParticles (const bool do_back_transformed_particles);
249  void SetDoBackTransformedParticles (std::string species_name, const bool do_back_transformed_particles);
250 
252  bool const onMainGrid = true;
253  auto const & v = m_deposit_on_main_grid;
254  return std::count( v.begin(), v.end(), onMainGrid );
255  }
256 
258  bool const fromMainGrid = true;
259  auto const & v = m_gather_from_main_grid;
260  return std::count( v.begin(), v.end(), fromMainGrid );
261  }
262 
263  // Inject particles during the simulation (for particles entering the
264  // simulation domain after some iterations, due to flowing plasma and/or
265  // moving window).
266  void ContinuousInjection(const amrex::RealBox& injection_box) const;
267  // Update injection position for continuously-injected species.
268  void UpdateContinuousInjectionPosition(amrex::Real dt) const;
269  int doContinuousInjection() const;
270 
271  // Inject particles from a surface during the simulation
272  void ContinuousFluxInjection(amrex::Real t, amrex::Real dt) const;
273 
274  std::vector<std::string> GetSpeciesNames() const { return species_names; }
275 
276  std::vector<std::string> GetLasersNames() const { return lasers_names; }
277 
278  std::vector<std::string> GetSpeciesAndLasersNames() const {
279  std::vector<std::string> tmp = species_names;
280  tmp.insert(tmp.end(), lasers_names.begin(), lasers_names.end());
281  return tmp;
282  }
283 
285 
286  void ScrapeParticles (const amrex::Vector<const amrex::MultiFab*>& distance_to_eb);
287 
288  std::string m_B_ext_particle_s = "none";
289  std::string m_E_ext_particle_s = "none";
290  // External fields added to particle fields.
293  // Parser for B_external on the particle
294  std::unique_ptr<amrex::Parser> m_Bx_particle_parser;
295  std::unique_ptr<amrex::Parser> m_By_particle_parser;
296  std::unique_ptr<amrex::Parser> m_Bz_particle_parser;
297  // Parser for E_external on the particle
298  std::unique_ptr<amrex::Parser> m_Ex_particle_parser;
299  std::unique_ptr<amrex::Parser> m_Ey_particle_parser;
300  std::unique_ptr<amrex::Parser> m_Ez_particle_parser;
301 
302  amrex::ParticleReal m_repeated_plasma_lens_period;
311 
312 #ifdef WARPX_QED
316  void doQedEvents (int lev,
317  const amrex::MultiFab& Ex,
318  const amrex::MultiFab& Ey,
319  const amrex::MultiFab& Ez,
320  const amrex::MultiFab& Bx,
321  const amrex::MultiFab& By,
322  const amrex::MultiFab& Bz);
323 #endif
324 
325  int getSpeciesID (std::string product_str) const;
326 
327 protected:
328 
329 #ifdef WARPX_QED
333  void doQedBreitWheeler (int lev,
334  const amrex::MultiFab& Ex,
335  const amrex::MultiFab& Ey,
336  const amrex::MultiFab& Ez,
337  const amrex::MultiFab& Bx,
338  const amrex::MultiFab& By,
339  const amrex::MultiFab& Bz);
340 
344  void doQedQuantumSync (int lev,
345  const amrex::MultiFab& Ex,
346  const amrex::MultiFab& Ey,
347  const amrex::MultiFab& Ez,
348  const amrex::MultiFab& Bx,
349  const amrex::MultiFab& By,
350  const amrex::MultiFab& Bz);
351 #endif
352 
353  // Particle container types
355 
356  std::vector<std::string> species_names;
357 
358  std::vector<std::string> lasers_names;
359 
360  //std::unique_ptr<CollisionHandler> collisionhandler;
361 
363  std::vector<bool> m_deposit_on_main_grid;
364  std::vector<bool> m_laser_deposit_on_main_grid;
365 
367  std::vector<bool> m_gather_from_main_grid;
368 
369  std::vector<PCTypes> species_types;
370 
371  template<typename ...Args>
373  Args const&... pc_dsts) const noexcept
374  {
375  amrex::MFItInfo info;
376 
377  MFItInfoCheckTiling(pc_src, pc_dsts...);
378 
379  if (pc_src.do_tiling && amrex::Gpu::notInLaunchRegion()) {
380  info.EnableTiling(pc_src.tile_size);
381  }
382 
383 #ifdef AMREX_USE_OMP
384  info.SetDynamic(true);
385 #endif
386 
387  return info;
388  }
389 
390 
391 #ifdef WARPX_QED
392  // The QED engines
393  std::shared_ptr<BreitWheelerEngine> m_shr_p_bw_engine;
394  std::shared_ptr<QuantumSynchrotronEngine> m_shr_p_qs_engine;
395  //_______________________________
396 
401  void InitQED ();
402 
403  //Variables to store how many species need a QED process
406  //________
407 
409  static_cast<amrex::ParticleReal>(
420 
425 
429  void InitQuantumSync ();
430 
434  void InitBreitWheeler ();
435 
441 
447 
468  amrex::Real m_qed_schwinger_xmin = std::numeric_limits<amrex::Real>::lowest();
469  amrex::Real m_qed_schwinger_xmax = std::numeric_limits<amrex::Real>::max();
470  amrex::Real m_qed_schwinger_ymin = std::numeric_limits<amrex::Real>::lowest();
471  amrex::Real m_qed_schwinger_ymax = std::numeric_limits<amrex::Real>::max();
472  amrex::Real m_qed_schwinger_zmin = std::numeric_limits<amrex::Real>::lowest();
473  amrex::Real m_qed_schwinger_zmax = std::numeric_limits<amrex::Real>::max();
474 
475 #endif
476 
477 private:
478 
479  // physical particles (+ laser)
481  // Temporary particle container, used e.g. for particle splitting.
482  std::unique_ptr<PhysicalParticleContainer> pc_tmp;
483 
484  void ReadParameters ();
485 
486  void mapSpeciesProduct ();
487 
489 
490  void MFItInfoCheckTiling(const WarpXParticleContainer& /*pc_src*/) const noexcept
491  {
492  return;
493  }
494 
495  template<typename First, typename ...Args>
497  First const& pc_dst, Args const&... others) const noexcept
498  {
499  if (pc_src.do_tiling && amrex::Gpu::notInLaunchRegion()) {
500  WARPX_ALWAYS_ASSERT_WITH_MESSAGE(pc_dst.do_tiling,
501  "For particle creation processes, either all or none of the "
502  "particle species must use tiling.");
503  }
504 
505  MFItInfoCheckTiling(pc_src, others...);
506  }
507 
514 
515 #ifdef WARPX_QED
521  void CheckQEDProductSpecies();
522 #endif
523 
524 
525 };
526 #endif /*WARPX_ParticleContainer_H_*/
#define WARPX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition: TextMsg.H:13
DtType
Definition: WarpXDtType.H:11
Definition: MultiParticleContainer.H:65
std::unique_ptr< amrex::Parser > m_Bz_particle_parser
Definition: MultiParticleContainer.H:296
WarpXParticleContainer * GetParticleContainerPtr(int ispecies) const
Definition: MultiParticleContainer.H:77
void InitQED()
Definition: MultiParticleContainer.cpp:979
std::unique_ptr< amrex::Parser > m_Ey_particle_parser
Definition: MultiParticleContainer.H:299
std::unique_ptr< amrex::MultiFab > GetChargeDensity(int lev, bool local=false)
Definition: MultiParticleContainer.cpp:602
amrex::Real m_qed_schwinger_zmax
Definition: MultiParticleContainer.H:473
void RedistributeLocal(const int num_ghost)
Definition: MultiParticleContainer.cpp:648
amrex::Real m_qed_schwinger_zmin
Definition: MultiParticleContainer.H:472
void InitMultiPhysicsModules()
Definition: MultiParticleContainer.cpp:452
std::string m_B_ext_particle_s
Definition: MultiParticleContainer.H:288
void QuantumSyncGenerateTable()
Definition: MultiParticleContainer.cpp:1139
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_strengths_E
Definition: MultiParticleContainer.H:305
void AllocData()
Definition: MultiParticleContainer.cpp:420
void SetParticleDistributionMap(int lev, amrex::DistributionMapping &new_dm)
Definition: MultiParticleContainer.cpp:712
void InitData()
Definition: MultiParticleContainer.cpp:429
amrex::Real m_qed_schwinger_ymin
Definition: MultiParticleContainer.H:470
std::vector< bool > m_laser_deposit_on_main_grid
Definition: MultiParticleContainer.H:364
void Increment(amrex::MultiFab &mf, int lev)
Definition: MultiParticleContainer.cpp:696
std::unique_ptr< WarpXParticleContainer > & GetUniqueContainer(int ispecies)
Definition: MultiParticleContainer.H:83
void mapSpeciesProduct()
Definition: MultiParticleContainer.cpp:778
std::unique_ptr< amrex::MultiFab > GetZeroChargeDensity(const int lev)
This returns a MultiFAB filled with zeros. It is used to return the charge density when there is no p...
Definition: MultiParticleContainer.cpp:515
void ContinuousInjection(const amrex::RealBox &injection_box) const
Definition: MultiParticleContainer.cpp:725
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_lengths
Definition: MultiParticleContainer.H:308
int NSpeciesBreitWheeler() const
Definition: MultiParticleContainer.H:424
std::array< amrex::ParticleReal, 3 > meanParticleVelocity(int ispecies)
Definition: MultiParticleContainer.H:87
void ReadHeader(std::istream &is)
Definition: ParticleIO.cpp:217
static constexpr auto m_default_quantum_sync_photon_creation_energy_threshold
Definition: MultiParticleContainer.H:408
void ApplyBoundaryConditions()
Definition: MultiParticleContainer.cpp:656
std::vector< std::string > lasers_names
Definition: MultiParticleContainer.H:358
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_starts
Definition: MultiParticleContainer.H:307
bool m_do_qed_schwinger
Definition: MultiParticleContainer.H:449
std::vector< std::string > GetLasersNames() const
Definition: MultiParticleContainer.H:276
std::unique_ptr< amrex::Parser > m_By_particle_parser
Definition: MultiParticleContainer.H:295
void PushX(amrex::Real dt)
Definition: MultiParticleContainer.cpp:497
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_lengths
Definition: MultiParticleContainer.H:304
std::string m_qed_schwinger_ele_product_name
Definition: MultiParticleContainer.H:451
int m_nspecies_quantum_sync
Definition: MultiParticleContainer.H:404
int m_qed_schwinger_ele_product
Definition: MultiParticleContainer.H:455
int nSpecies() const
Definition: MultiParticleContainer.H:237
int m_nspecies_breit_wheeler
Definition: MultiParticleContainer.H:405
void MFItInfoCheckTiling(const WarpXParticleContainer &) const noexcept
Definition: MultiParticleContainer.H:490
void InitQuantumSync()
Definition: MultiParticleContainer.cpp:1008
std::shared_ptr< BreitWheelerEngine > m_shr_p_bw_engine
Definition: MultiParticleContainer.H:393
amrex::Vector< amrex::ParticleReal > m_E_external_particle
Definition: MultiParticleContainer.H:292
void Redistribute()
Definition: MultiParticleContainer.cpp:632
std::vector< PCTypes > species_types
Definition: MultiParticleContainer.H:369
void PostRestart()
Definition: MultiParticleContainer.cpp:441
int NSpeciesQuantumSync() const
Definition: MultiParticleContainer.H:419
amrex::ParticleReal m_quantum_sync_photon_creation_energy_threshold
Definition: MultiParticleContainer.H:413
WarpXParticleContainer & GetParticleContainerFromName(const std::string &name) const
Definition: MultiParticleContainer.cpp:409
std::vector< std::string > species_names
Definition: MultiParticleContainer.H:356
amrex::Box ComputeSchwingerGlobalBox() const
Definition: MultiParticleContainer.cpp:1417
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_strengths_B
Definition: MultiParticleContainer.H:310
amrex::Real m_qed_schwinger_ymax
Definition: MultiParticleContainer.H:471
void doQEDSchwinger()
Definition: MultiParticleContainer.cpp:1312
std::string m_E_ext_particle_s
Definition: MultiParticleContainer.H:289
amrex::Real m_qed_schwinger_xmin
Definition: MultiParticleContainer.H:468
void UpdateContinuousInjectionPosition(amrex::Real dt) const
Definition: MultiParticleContainer.cpp:740
int nSpeciesGatherFromMainGrid() const
Definition: MultiParticleContainer.H:257
void defineAllParticleTiles()
Definition: MultiParticleContainer.cpp:640
void BreitWheelerGenerateTable()
Definition: MultiParticleContainer.cpp:1228
std::unique_ptr< PhysicalParticleContainer > pc_tmp
Definition: MultiParticleContainer.H:482
void ReadParameters()
Definition: MultiParticleContainer.cpp:128
int doContinuousInjection() const
Definition: MultiParticleContainer.cpp:750
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_strengths_B
Definition: MultiParticleContainer.H:306
int nSpeciesDepositOnMainGrid() const
Definition: MultiParticleContainer.H:251
amrex::ParticleReal m_repeated_plasma_lens_period
Definition: MultiParticleContainer.H:302
amrex::Real m_qed_schwinger_xmax
Definition: MultiParticleContainer.H:469
void PushP(int lev, amrex::Real dt, const amrex::MultiFab &Ex, const amrex::MultiFab &Ey, const amrex::MultiFab &Ez, const amrex::MultiFab &Bx, const amrex::MultiFab &By, const amrex::MultiFab &Bz)
Definition: MultiParticleContainer.cpp:505
amrex::Vector< amrex::ParticleReal > m_B_external_particle
Definition: MultiParticleContainer.H:291
void MFItInfoCheckTiling(const WarpXParticleContainer &pc_src, First const &pc_dst, Args const &... others) const noexcept
Definition: MultiParticleContainer.H:496
int getSpeciesID(std::string product_str) const
Definition: MultiParticleContainer.cpp:823
std::vector< std::string > GetSpeciesAndLasersNames() const
Definition: MultiParticleContainer.H:278
~MultiParticleContainer()
Definition: MultiParticleContainer.H:71
amrex::Vector< amrex::Long > NumberOfParticlesInGrid(int lev) const
Definition: MultiParticleContainer.cpp:673
std::vector< std::string > GetSpeciesNames() const
Definition: MultiParticleContainer.H:274
void InitBreitWheeler()
Definition: MultiParticleContainer.cpp:1080
void doFieldIonization(int lev, const amrex::MultiFab &Ex, const amrex::MultiFab &Ey, const amrex::MultiFab &Ez, const amrex::MultiFab &Bx, const amrex::MultiFab &By, const amrex::MultiFab &Bz)
Definition: MultiParticleContainer.cpp:871
bool m_do_back_transformed_particles
Definition: MultiParticleContainer.H:488
void SortParticlesByBin(amrex::IntVect bin_size)
Definition: MultiParticleContainer.cpp:620
amrex::Vector< std::unique_ptr< WarpXParticleContainer > > allcontainers
Definition: MultiParticleContainer.H:480
void SetParticleBoxArray(int lev, amrex::BoxArray &new_ba)
Definition: MultiParticleContainer.cpp:704
void doQedQuantumSync(int lev, const amrex::MultiFab &Ex, const amrex::MultiFab &Ey, const amrex::MultiFab &Ez, const amrex::MultiFab &Bx, const amrex::MultiFab &By, const amrex::MultiFab &Bz)
Performs QED photon emission for the species for which it is enabled.
Definition: MultiParticleContainer.cpp:1574
void CheckIonizationProductSpecies()
Definition: MultiParticleContainer.cpp:956
std::unique_ptr< amrex::Parser > m_Ex_particle_parser
Definition: MultiParticleContainer.H:298
amrex::Vector< amrex::Long > GetZeroParticlesInGrid(const int lev) const
This returns a vector filled with zeros whose size is the number of boxes in the simulation boxarray....
Definition: MultiParticleContainer.cpp:664
void doResampling(const int timestep)
This function loops over all species and performs resampling if appropriate.
Definition: MultiParticleContainer.cpp:945
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_starts
Definition: MultiParticleContainer.H:303
int m_qed_schwinger_threshold_poisson_gaussian
Definition: MultiParticleContainer.H:464
void SetDoBackTransformedParticles(const bool do_back_transformed_particles)
Definition: MultiParticleContainer.cpp:847
std::vector< bool > m_gather_from_main_grid
instead of gathering fields from the finest patch level, gather from the coarsest
Definition: MultiParticleContainer.H:367
void doQedEvents(int lev, const amrex::MultiFab &Ex, const amrex::MultiFab &Ey, const amrex::MultiFab &Ez, const amrex::MultiFab &Bx, const amrex::MultiFab &By, const amrex::MultiFab &Bz)
Performs QED events (Breit-Wheeler process and photon emission)
Definition: MultiParticleContainer.cpp:1479
void CheckQEDProductSpecies()
Definition: MultiParticleContainer.cpp:1653
amrex::Real m_qed_schwinger_y_size
Definition: MultiParticleContainer.H:459
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_strengths_E
Definition: MultiParticleContainer.H:309
void ScrapeParticles(const amrex::Vector< const amrex::MultiFab * > &distance_to_eb)
Definition: MultiParticleContainer.cpp:967
PhysicalParticleContainer & GetPCtmp()
Definition: MultiParticleContainer.H:284
void DepositCurrent(amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > &J, const amrex::Real dt, const amrex::Real relative_time)
Deposit current density.
Definition: MultiParticleContainer.cpp:537
void WriteHeader(std::ostream &os) const
Definition: ParticleIO.cpp:228
std::string m_qed_schwinger_pos_product_name
Definition: MultiParticleContainer.H:453
void Restart(const std::string &dir)
Definition: ParticleIO.cpp:121
WarpXParticleContainer & GetParticleContainer(int ispecies) const
Definition: MultiParticleContainer.H:74
void ContinuousFluxInjection(amrex::Real t, amrex::Real dt) const
Definition: MultiParticleContainer.cpp:766
void DepositCharge(amrex::Vector< std::unique_ptr< amrex::MultiFab > > &rho, const amrex::Real relative_time)
Deposit charge density.
Definition: MultiParticleContainer.cpp:565
int m_qed_schwinger_pos_product
Definition: MultiParticleContainer.H:457
amrex::MFItInfo getMFItInfo(const WarpXParticleContainer &pc_src, Args const &... pc_dsts) const noexcept
Definition: MultiParticleContainer.H:372
MultiParticleContainer(amrex::AmrCore *amr_core)
Definition: MultiParticleContainer.cpp:92
std::unique_ptr< amrex::Parser > m_Ez_particle_parser
Definition: MultiParticleContainer.H:300
std::shared_ptr< QuantumSynchrotronEngine > m_shr_p_qs_engine
Definition: MultiParticleContainer.H:394
std::vector< bool > m_deposit_on_main_grid
instead of depositing (current, charge) on the finest patch level, deposit to the coarsest grid
Definition: MultiParticleContainer.H:363
PCTypes
Definition: MultiParticleContainer.H:354
void Evolve(int lev, const amrex::MultiFab &Ex, const amrex::MultiFab &Ey, const amrex::MultiFab &Ez, const amrex::MultiFab &Bx, const amrex::MultiFab &By, const amrex::MultiFab &Bz, amrex::MultiFab &jx, amrex::MultiFab &jy, amrex::MultiFab &jz, amrex::MultiFab *cjx, amrex::MultiFab *cjy, amrex::MultiFab *cjz, amrex::MultiFab *rho, amrex::MultiFab *crho, const amrex::MultiFab *cEx, const amrex::MultiFab *cEy, const amrex::MultiFab *cEz, const amrex::MultiFab *cBx, const amrex::MultiFab *cBy, const amrex::MultiFab *cBz, amrex::Real t, amrex::Real dt, DtType a_dt_type=DtType::Full, bool skip_deposition=false)
Definition: MultiParticleContainer.cpp:470
std::unique_ptr< amrex::Parser > m_Bx_particle_parser
Definition: MultiParticleContainer.H:294
void doQedBreitWheeler(int lev, const amrex::MultiFab &Ex, const amrex::MultiFab &Ey, const amrex::MultiFab &Ez, const amrex::MultiFab &Bx, const amrex::MultiFab &By, const amrex::MultiFab &Bz)
Performs Breit-Wheeler process for the species for which it is enabled.
Definition: MultiParticleContainer.cpp:1493
Definition: PhysicalParticleContainer.H:46
Definition: WarpXParticleContainer.H:104
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
bool notInLaunchRegion() noexcept
int count
Definition: run_alltests.py:322
string name
Definition: stencil.py:452
int dt
Definition: stencil.py:440
MFItInfo & SetDynamic(bool f) noexcept
MFItInfo & EnableTiling(const IntVect &ts=FabArrayBase::mfiter_tile_size) noexcept