WarpX
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"
17 #include "Evolve/WarpXPushType.H"
19 #ifdef WARPX_QED
22 #endif
24 #include "Utils/TextMsg.H"
25 #include "Utils/WarpXConst.H"
26 #include "WarpXParticleContainer.H"
27 #include "ParticleBoundaries.H"
28 
29 #include <AMReX_BLassert.H>
30 #include <AMReX_Box.H>
31 #include <AMReX_Config.H>
32 #include <AMReX_GpuControl.H>
33 #include <AMReX_INT.H>
34 #include <AMReX_MFIter.H>
35 #include <AMReX_REAL.H>
36 #include <AMReX_RealBox.H>
37 #include <AMReX_Vector.H>
38 
39 #include <AMReX_BaseFwd.H>
40 #include <AMReX_AmrCoreFwd.H>
41 
42 #include <algorithm>
43 #include <array>
44 #include <iosfwd>
45 #include <iterator>
46 #include <limits>
47 #include <memory>
48 #include <string>
49 #include <vector>
50 
66 {
67 
68 public:
69 
71 
73 
78 
79  [[nodiscard]] WarpXParticleContainer&
81 
82  [[nodiscard]] WarpXParticleContainer*
83  GetParticleContainerPtr (int index) const {return allcontainers[index].get();}
84 
85  [[nodiscard]] WarpXParticleContainer&
86  GetParticleContainerFromName (const std::string& name) const;
87 
88  std::array<amrex::ParticleReal, 3> meanParticleVelocity(int index) {
89  return allcontainers[index]->meanParticleVelocity();
90  }
91 
92  void AllocData ();
93 
94  void InitData ();
95 
97 
103  void Evolve (int lev,
104  const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
105  const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz,
108  amrex::MultiFab* rho, amrex::MultiFab* crho,
109  const amrex::MultiFab* cEx, const amrex::MultiFab* cEy, const amrex::MultiFab* cEz,
110  const amrex::MultiFab* cBx, const amrex::MultiFab* cBy, const amrex::MultiFab* cBz,
111  amrex::Real t, amrex::Real dt, DtType a_dt_type=DtType::Full, bool skip_deposition=false,
112  PushType push_type=PushType::Explicit);
113 
119  void PushX (amrex::Real dt);
120 
127  void PushP (int lev, amrex::Real dt,
128  const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
129  const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz);
130 
137  std::unique_ptr<amrex::MultiFab> GetZeroChargeDensity(int lev);
138 
148  void
149  DepositCharge (amrex::Vector<std::unique_ptr<amrex::MultiFab> >& rho,
150  amrex::Real relative_time);
151 
163  void
164  DepositCurrent (amrex::Vector<std::array< std::unique_ptr<amrex::MultiFab>, 3 > >& J,
165  amrex::Real dt, amrex::Real relative_time);
166 
170 
171  std::unique_ptr<amrex::MultiFab> GetChargeDensity(int lev, bool local = false);
172 
173  void doFieldIonization (int lev,
174  const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
175  const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz);
176 
177  void doCollisions (amrex::Real cur_time, amrex::Real dt);
178 
184  void doResampling (int timestep, bool verbose);
185 
186 #ifdef WARPX_QED
194  void doQEDSchwinger ();
195 
200  [[nodiscard]] amrex::Box ComputeSchwingerGlobalBox () const;
201 #endif
202 
203  void Restart (const std::string& dir);
204 
205  void PostRestart ();
206 
207  void ReadHeader (std::istream& is);
208 
209  void WriteHeader (std::ostream& os) const;
210 
211  void SortParticlesByBin (amrex::IntVect bin_size);
212 
213  void Redistribute ();
214 
215  void defineAllParticleTiles ();
216 
217  void RedistributeLocal (int num_ghost);
218 
221  void ApplyBoundaryConditions ();
222 
230  [[nodiscard]] amrex::Vector<amrex::Long> GetZeroParticlesInGrid(int lev) const;
231 
232  [[nodiscard]] amrex::Vector<amrex::Long> NumberOfParticlesInGrid(int lev) const;
233 
234  void Increment (amrex::MultiFab& mf, int lev);
235 
236  void SetParticleBoxArray (int lev, amrex::BoxArray& new_ba);
238 
239  [[nodiscard]] int nSpecies () const {return static_cast<int>(species_names.size());}
240  [[nodiscard]] int nLasers () const {return static_cast<int>(lasers_names.size());}
241  [[nodiscard]] int nContainers () const {return static_cast<int>(allcontainers.size());}
242 
247  void SetDoBackTransformedParticles (bool do_back_transformed_particles);
253  void SetDoBackTransformedParticles (std::string species_name, bool do_back_transformed_particles);
254 
255  [[nodiscard]] int nSpeciesDepositOnMainGrid () const
256  {
257  bool const onMainGrid = true;
258  auto const & v = m_deposit_on_main_grid;
259  return static_cast<int>(std::count( v.begin(), v.end(), onMainGrid ));
260  }
261 
262  [[nodiscard]] int nSpeciesGatherFromMainGrid() const
263  {
264  bool const fromMainGrid = true;
265  auto const & v = m_gather_from_main_grid;
266  return static_cast<int>(std::count( v.begin(), v.end(), fromMainGrid ));
267  }
268 
269  // Inject particles during the simulation (for particles entering the
270  // simulation domain after some iterations, due to flowing plasma and/or
271  // moving window).
272  void ContinuousInjection(const amrex::RealBox& injection_box) const;
273 
278  void UpdateAntennaPosition(amrex::Real dt) const;
279 
280  [[nodiscard]] int doContinuousInjection() const;
281 
282  // Inject particles from a surface during the simulation
283  void ContinuousFluxInjection(amrex::Real t, amrex::Real dt) const;
284 
285  [[nodiscard]] std::vector<std::string> GetSpeciesNames() const { return species_names; }
286 
287  [[nodiscard]] std::vector<std::string> GetLasersNames() const { return lasers_names; }
288 
289  [[nodiscard]] std::vector<std::string> GetSpeciesAndLasersNames() const
290  {
291  std::vector<std::string> tmp = species_names;
292  tmp.insert(tmp.end(), lasers_names.begin(), lasers_names.end());
293  return tmp;
294  }
295 
297 
298  void ScrapeParticles (const amrex::Vector<const amrex::MultiFab*>& distance_to_eb);
299 
300  std::string m_B_ext_particle_s = "none";
301  std::string m_E_ext_particle_s = "none";
302  // Parser for B_external on the particle
303  std::unique_ptr<amrex::Parser> m_Bx_particle_parser;
304  std::unique_ptr<amrex::Parser> m_By_particle_parser;
305  std::unique_ptr<amrex::Parser> m_Bz_particle_parser;
306  // Parser for E_external on the particle
307  std::unique_ptr<amrex::Parser> m_Ex_particle_parser;
308  std::unique_ptr<amrex::Parser> m_Ey_particle_parser;
309  std::unique_ptr<amrex::Parser> m_Ez_particle_parser;
310 
311  amrex::ParticleReal m_repeated_plasma_lens_period;
320 
321 #ifdef WARPX_QED
325  void doQedEvents (int lev,
326  const amrex::MultiFab& Ex,
327  const amrex::MultiFab& Ey,
328  const amrex::MultiFab& Ez,
329  const amrex::MultiFab& Bx,
330  const amrex::MultiFab& By,
331  const amrex::MultiFab& Bz);
332 #endif
333 
334  [[nodiscard]] int getSpeciesID (std::string product_str) const;
335 
338 
339 protected:
340 
341 #ifdef WARPX_QED
345  void doQedBreitWheeler (int lev,
346  const amrex::MultiFab& Ex,
347  const amrex::MultiFab& Ey,
348  const amrex::MultiFab& Ez,
349  const amrex::MultiFab& Bx,
350  const amrex::MultiFab& By,
351  const amrex::MultiFab& Bz);
352 
356  void doQedQuantumSync (int lev,
357  const amrex::MultiFab& Ex,
358  const amrex::MultiFab& Ey,
359  const amrex::MultiFab& Ez,
360  const amrex::MultiFab& Bx,
361  const amrex::MultiFab& By,
362  const amrex::MultiFab& Bz);
363 #endif
364 
365  // Particle container types
367 
368  std::vector<std::string> species_names;
369 
370  std::vector<std::string> lasers_names;
371 
372  std::unique_ptr<CollisionHandler> collisionhandler;
373 
375  std::vector<bool> m_deposit_on_main_grid;
376  std::vector<bool> m_laser_deposit_on_main_grid;
377 
379  std::vector<bool> m_gather_from_main_grid;
380 
381  std::vector<PCTypes> species_types;
382 
383  template<typename ...Args>
384  [[nodiscard]]
386  Args const&... pc_dsts) const noexcept
387  {
388  amrex::MFItInfo info;
389 
390  MFItInfoCheckTiling(pc_src, pc_dsts...);
391 
394  }
395 
396 #ifdef AMREX_USE_OMP
397  info.SetDynamic(true);
398 #endif
399 
400  return info;
401  }
402 
403 
404 #ifdef WARPX_QED
405  // The QED engines
406  std::shared_ptr<BreitWheelerEngine> m_shr_p_bw_engine;
407  std::shared_ptr<QuantumSynchrotronEngine> m_shr_p_qs_engine;
408  //_______________________________
409 
414  void InitQED ();
415 
416  //Variables to store how many species need a QED process
419  //________
420 
422  static_cast<amrex::ParticleReal>(
432  [[nodiscard]] int NSpeciesQuantumSync() const { return m_nspecies_quantum_sync;}
433 
437  [[nodiscard]] int NSpeciesBreitWheeler() const { return m_nspecies_breit_wheeler;}
438 
442  void InitQuantumSync ();
443 
447  void InitBreitWheeler ();
448 
454 
460 
462  bool m_do_qed_schwinger = false;
481  amrex::Real m_qed_schwinger_xmin = std::numeric_limits<amrex::Real>::lowest();
482  amrex::Real m_qed_schwinger_xmax = std::numeric_limits<amrex::Real>::max();
483  amrex::Real m_qed_schwinger_ymin = std::numeric_limits<amrex::Real>::lowest();
484  amrex::Real m_qed_schwinger_ymax = std::numeric_limits<amrex::Real>::max();
485  amrex::Real m_qed_schwinger_zmin = std::numeric_limits<amrex::Real>::lowest();
486  amrex::Real m_qed_schwinger_zmax = std::numeric_limits<amrex::Real>::max();
487 
488 #endif
489 
490 private:
491 
492  // physical particles (+ laser)
494  // Temporary particle container, used e.g. for particle splitting.
495  std::unique_ptr<PhysicalParticleContainer> pc_tmp;
496 
497  void ReadParameters ();
498 
499  void mapSpeciesProduct ();
500 
502 
503  void MFItInfoCheckTiling(const WarpXParticleContainer& /*pc_src*/) const noexcept
504  {}
505 
506  template<typename First, typename ...Args>
508  First const& pc_dst, Args const&... others) const noexcept
509  {
511  WARPX_ALWAYS_ASSERT_WITH_MESSAGE(pc_dst.do_tiling,
512  "For particle creation processes, either all or none of the "
513  "particle species must use tiling.");
514  }
515 
516  MFItInfoCheckTiling(pc_src, others...);
517  }
518 
525 
526 #ifdef WARPX_QED
532  void CheckQEDProductSpecies();
533 #endif
534 
535 
536 };
537 #endif /*WARPX_ParticleContainer_H_*/
#define WARPX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition: TextMsg.H:13
DtType
Definition: WarpXDtType.H:11
PushType
Definition: WarpXPushType.H:12
Definition: MultiParticleContainer.H:66
void RedistributeLocal(int num_ghost)
Definition: MultiParticleContainer.cpp:636
std::unique_ptr< amrex::Parser > m_Bz_particle_parser
Definition: MultiParticleContainer.H:305
void InitQED()
Definition: MultiParticleContainer.cpp:960
std::unique_ptr< amrex::Parser > m_Ey_particle_parser
Definition: MultiParticleContainer.H:308
std::unique_ptr< amrex::MultiFab > GetChargeDensity(int lev, bool local=false)
Definition: MultiParticleContainer.cpp:586
int nLasers() const
Definition: MultiParticleContainer.H:240
amrex::Real m_qed_schwinger_zmax
Definition: MultiParticleContainer.H:486
amrex::Real m_qed_schwinger_zmin
Definition: MultiParticleContainer.H:485
void InitMultiPhysicsModules()
Definition: MultiParticleContainer.cpp:432
std::string m_B_ext_particle_s
Definition: MultiParticleContainer.H:300
amrex::Vector< std::unique_ptr< WarpXParticleContainer > >::iterator begin()
Definition: MultiParticleContainer.H:336
void QuantumSyncGenerateTable()
Definition: MultiParticleContainer.cpp:1123
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_strengths_E
Definition: MultiParticleContainer.H:314
void AllocData()
Definition: MultiParticleContainer.cpp:400
void SetParticleDistributionMap(int lev, amrex::DistributionMapping &new_dm)
Definition: MultiParticleContainer.cpp:700
void InitData()
Definition: MultiParticleContainer.cpp:409
amrex::Real m_qed_schwinger_ymin
Definition: MultiParticleContainer.H:483
std::vector< bool > m_laser_deposit_on_main_grid
Definition: MultiParticleContainer.H:376
void Increment(amrex::MultiFab &mf, int lev)
Definition: MultiParticleContainer.cpp:684
void mapSpeciesProduct()
Definition: MultiParticleContainer.cpp:761
void ContinuousInjection(const amrex::RealBox &injection_box) const
Definition: MultiParticleContainer.cpp:713
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_lengths
Definition: MultiParticleContainer.H:317
MultiParticleContainer(MultiParticleContainer const &)=delete
int NSpeciesBreitWheeler() const
Definition: MultiParticleContainer.H:437
void ReadHeader(std::istream &is)
Definition: ParticleIO.cpp:218
static constexpr auto m_default_quantum_sync_photon_creation_energy_threshold
Definition: MultiParticleContainer.H:421
void ApplyBoundaryConditions()
Definition: MultiParticleContainer.cpp:644
std::vector< std::string > lasers_names
Definition: MultiParticleContainer.H:370
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_starts
Definition: MultiParticleContainer.H:316
bool m_do_qed_schwinger
Definition: MultiParticleContainer.H:462
std::vector< std::string > GetLasersNames() const
Definition: MultiParticleContainer.H:287
std::unique_ptr< amrex::Parser > m_By_particle_parser
Definition: MultiParticleContainer.H:304
void PushX(amrex::Real dt)
Definition: MultiParticleContainer.cpp:478
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_lengths
Definition: MultiParticleContainer.H:313
std::string m_qed_schwinger_ele_product_name
Definition: MultiParticleContainer.H:464
int m_nspecies_quantum_sync
Definition: MultiParticleContainer.H:417
int m_qed_schwinger_ele_product
Definition: MultiParticleContainer.H:468
int nSpecies() const
Definition: MultiParticleContainer.H:239
int m_nspecies_breit_wheeler
Definition: MultiParticleContainer.H:418
void MFItInfoCheckTiling(const WarpXParticleContainer &) const noexcept
Definition: MultiParticleContainer.H:503
void InitQuantumSync()
Definition: MultiParticleContainer.cpp:991
void DepositCurrent(amrex::Vector< std::array< std::unique_ptr< amrex::MultiFab >, 3 > > &J, amrex::Real dt, amrex::Real relative_time)
Deposit current density.
Definition: MultiParticleContainer.cpp:521
std::shared_ptr< BreitWheelerEngine > m_shr_p_bw_engine
Definition: MultiParticleContainer.H:406
void Redistribute()
Definition: MultiParticleContainer.cpp:620
std::vector< PCTypes > species_types
Definition: MultiParticleContainer.H:381
void PostRestart()
Definition: MultiParticleContainer.cpp:421
int NSpeciesQuantumSync() const
Definition: MultiParticleContainer.H:432
amrex::ParticleReal m_quantum_sync_photon_creation_energy_threshold
Definition: MultiParticleContainer.H:426
WarpXParticleContainer & GetParticleContainerFromName(const std::string &name) const
Definition: MultiParticleContainer.cpp:389
MultiParticleContainer & operator=(MultiParticleContainer &&)=default
std::vector< std::string > species_names
Definition: MultiParticleContainer.H:368
amrex::Box ComputeSchwingerGlobalBox() const
Definition: MultiParticleContainer.cpp:1403
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_strengths_B
Definition: MultiParticleContainer.H:319
amrex::Real m_qed_schwinger_ymax
Definition: MultiParticleContainer.H:484
void doQEDSchwinger()
Definition: MultiParticleContainer.cpp:1298
std::string m_E_ext_particle_s
Definition: MultiParticleContainer.H:301
amrex::Real m_qed_schwinger_xmin
Definition: MultiParticleContainer.H:481
int nSpeciesGatherFromMainGrid() const
Definition: MultiParticleContainer.H:262
void defineAllParticleTiles()
Definition: MultiParticleContainer.cpp:628
void BreitWheelerGenerateTable()
Definition: MultiParticleContainer.cpp:1213
std::unique_ptr< PhysicalParticleContainer > pc_tmp
Definition: MultiParticleContainer.H:495
std::unique_ptr< amrex::MultiFab > GetZeroChargeDensity(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:496
void ReadParameters()
Definition: MultiParticleContainer.cpp:129
int doContinuousInjection() const
Definition: MultiParticleContainer.cpp:733
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_strengths_B
Definition: MultiParticleContainer.H:315
int nSpeciesDepositOnMainGrid() const
Definition: MultiParticleContainer.H:255
void DepositCharge(amrex::Vector< std::unique_ptr< amrex::MultiFab > > &rho, amrex::Real relative_time)
Deposit charge density.
Definition: MultiParticleContainer.cpp:549
std::unique_ptr< CollisionHandler > collisionhandler
Definition: MultiParticleContainer.H:372
amrex::ParticleReal m_repeated_plasma_lens_period
Definition: MultiParticleContainer.H:311
amrex::Real m_qed_schwinger_xmax
Definition: MultiParticleContainer.H:482
WarpXParticleContainer * GetParticleContainerPtr(int index) const
Definition: MultiParticleContainer.H:83
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:486
void MFItInfoCheckTiling(const WarpXParticleContainer &pc_src, First const &pc_dst, Args const &... others) const noexcept
Definition: MultiParticleContainer.H:507
int getSpeciesID(std::string product_str) const
Definition: MultiParticleContainer.cpp:806
std::vector< std::string > GetSpeciesAndLasersNames() const
Definition: MultiParticleContainer.H:289
amrex::Vector< amrex::Long > NumberOfParticlesInGrid(int lev) const
Definition: MultiParticleContainer.cpp:661
std::vector< std::string > GetSpeciesNames() const
Definition: MultiParticleContainer.H:285
void InitBreitWheeler()
Definition: MultiParticleContainer.cpp:1063
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:854
amrex::Vector< amrex::Long > GetZeroParticlesInGrid(int lev) const
This returns a vector filled with zeros whose size is the number of boxes in the simulation boxarray....
Definition: MultiParticleContainer.cpp:652
bool m_do_back_transformed_particles
Definition: MultiParticleContainer.H:501
void SortParticlesByBin(amrex::IntVect bin_size)
Definition: MultiParticleContainer.cpp:608
amrex::Vector< std::unique_ptr< WarpXParticleContainer > > allcontainers
Definition: MultiParticleContainer.H:493
void SetParticleBoxArray(int lev, amrex::BoxArray &new_ba)
Definition: MultiParticleContainer.cpp:692
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:1562
void UpdateAntennaPosition(amrex::Real dt) const
Update antenna position for continuous injection of lasers in a boosted frame. Empty function for con...
Definition: MultiParticleContainer.cpp:723
MultiParticleContainer(MultiParticleContainer &&)=default
void doResampling(int timestep, bool verbose)
This function loops over all species and performs resampling if appropriate.
Definition: MultiParticleContainer.cpp:926
void CheckIonizationProductSpecies()
Definition: MultiParticleContainer.cpp:937
std::unique_ptr< amrex::Parser > m_Ex_particle_parser
Definition: MultiParticleContainer.H:307
std::array< amrex::ParticleReal, 3 > meanParticleVelocity(int index)
Definition: MultiParticleContainer.H:88
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_starts
Definition: MultiParticleContainer.H:312
int m_qed_schwinger_threshold_poisson_gaussian
Definition: MultiParticleContainer.H:477
std::vector< bool > m_gather_from_main_grid
instead of gathering fields from the finest patch level, gather from the coarsest
Definition: MultiParticleContainer.H:379
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:1465
void CheckQEDProductSpecies()
Definition: MultiParticleContainer.cpp:1643
MultiParticleContainer & operator=(MultiParticleContainer const &)=delete
amrex::Real m_qed_schwinger_y_size
Definition: MultiParticleContainer.H:472
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_strengths_E
Definition: MultiParticleContainer.H:318
void ScrapeParticles(const amrex::Vector< const amrex::MultiFab * > &distance_to_eb)
Definition: MultiParticleContainer.cpp:948
PhysicalParticleContainer & GetPCtmp()
Definition: MultiParticleContainer.H:296
int nContainers() const
Definition: MultiParticleContainer.H:241
void WriteHeader(std::ostream &os) const
Definition: ParticleIO.cpp:229
std::string m_qed_schwinger_pos_product_name
Definition: MultiParticleContainer.H:466
void Restart(const std::string &dir)
Definition: ParticleIO.cpp:122
void ContinuousFluxInjection(amrex::Real t, amrex::Real dt) const
Definition: MultiParticleContainer.cpp:749
int m_qed_schwinger_pos_product
Definition: MultiParticleContainer.H:470
amrex::MFItInfo getMFItInfo(const WarpXParticleContainer &pc_src, Args const &... pc_dsts) const noexcept
Definition: MultiParticleContainer.H:385
MultiParticleContainer(amrex::AmrCore *amr_core)
Definition: MultiParticleContainer.cpp:93
WarpXParticleContainer & GetParticleContainer(int index) const
Definition: MultiParticleContainer.H:80
void doCollisions(amrex::Real cur_time, amrex::Real dt)
Definition: MultiParticleContainer.cpp:920
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, PushType push_type=PushType::Explicit)
Definition: MultiParticleContainer.cpp:450
std::unique_ptr< amrex::Parser > m_Ez_particle_parser
Definition: MultiParticleContainer.H:309
std::shared_ptr< QuantumSynchrotronEngine > m_shr_p_qs_engine
Definition: MultiParticleContainer.H:407
~MultiParticleContainer()=default
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:375
void SetDoBackTransformedParticles(bool do_back_transformed_particles)
Definition: MultiParticleContainer.cpp:830
PCTypes
Definition: MultiParticleContainer.H:366
amrex::Vector< std::unique_ptr< WarpXParticleContainer > >::iterator end()
Definition: MultiParticleContainer.H:337
std::unique_ptr< amrex::Parser > m_Bx_particle_parser
Definition: MultiParticleContainer.H:303
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:1479
Definition: PhysicalParticleContainer.H:47
Definition: WarpXParticleContainer.H:110
static AMREX_EXPORT bool do_tiling
static AMREX_EXPORT IntVect tile_size
Long size() const noexcept
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
name
Definition: run_automated.py:229
index
Definition: run_automated.py:328
float dt
Definition: stencil.py:442
MFItInfo & EnableTiling(const IntVect &ts=FabArrayBase::mfiter_tile_size) noexcept
MFItInfo & SetDynamic(bool f) noexcept