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"
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 
77 
80 
82  GetParticleContainerPtr (int index) const {return allcontainers[index].get();}
83 
85  GetParticleContainerFromName (const std::string& name) const;
86 
87  std::array<amrex::ParticleReal, 3> meanParticleVelocity(int index) {
88  return allcontainers[index]->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(int lev);
136 
146  void
147  DepositCharge (amrex::Vector<std::unique_ptr<amrex::MultiFab> >& rho,
148  amrex::Real relative_time);
149 
161  void
162  DepositCurrent (amrex::Vector<std::array< std::unique_ptr<amrex::MultiFab>, 3 > >& J,
163  amrex::Real dt, 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 (int timestep, bool verbose);
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 (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 static_cast<int>(species_names.size());}
238  int nLasers () const {return static_cast<int>(lasers_names.size());}
239  int nContainers () const {return static_cast<int>(allcontainers.size());}
240 
245  void SetDoBackTransformedParticles (bool do_back_transformed_particles);
251  void SetDoBackTransformedParticles (std::string species_name, bool do_back_transformed_particles);
252 
254  bool const onMainGrid = true;
255  auto const & v = m_deposit_on_main_grid;
256  return static_cast<int>(std::count( v.begin(), v.end(), onMainGrid ));
257  }
258 
260  bool const fromMainGrid = true;
261  auto const & v = m_gather_from_main_grid;
262  return static_cast<int>(std::count( v.begin(), v.end(), fromMainGrid ));
263  }
264 
265  // Inject particles during the simulation (for particles entering the
266  // simulation domain after some iterations, due to flowing plasma and/or
267  // moving window).
268  void ContinuousInjection(const amrex::RealBox& injection_box) const;
269 
274  void UpdateAntennaPosition(amrex::Real dt) const;
275 
276  int doContinuousInjection() const;
277 
278  // Inject particles from a surface during the simulation
279  void ContinuousFluxInjection(amrex::Real t, amrex::Real dt) const;
280 
281  std::vector<std::string> GetSpeciesNames() const { return species_names; }
282 
283  std::vector<std::string> GetLasersNames() const { return lasers_names; }
284 
285  std::vector<std::string> GetSpeciesAndLasersNames() const {
286  std::vector<std::string> tmp = species_names;
287  tmp.insert(tmp.end(), lasers_names.begin(), lasers_names.end());
288  return tmp;
289  }
290 
292 
293  void ScrapeParticles (const amrex::Vector<const amrex::MultiFab*>& distance_to_eb);
294 
295  std::string m_B_ext_particle_s = "none";
296  std::string m_E_ext_particle_s = "none";
297  // External fields added to particle fields.
300  // Parser for B_external on the particle
301  std::unique_ptr<amrex::Parser> m_Bx_particle_parser;
302  std::unique_ptr<amrex::Parser> m_By_particle_parser;
303  std::unique_ptr<amrex::Parser> m_Bz_particle_parser;
304  // Parser for E_external on the particle
305  std::unique_ptr<amrex::Parser> m_Ex_particle_parser;
306  std::unique_ptr<amrex::Parser> m_Ey_particle_parser;
307  std::unique_ptr<amrex::Parser> m_Ez_particle_parser;
308 
309  amrex::ParticleReal m_repeated_plasma_lens_period;
318 
319 #ifdef WARPX_QED
323  void doQedEvents (int lev,
324  const amrex::MultiFab& Ex,
325  const amrex::MultiFab& Ey,
326  const amrex::MultiFab& Ez,
327  const amrex::MultiFab& Bx,
328  const amrex::MultiFab& By,
329  const amrex::MultiFab& Bz);
330 #endif
331 
332  int getSpeciesID (std::string product_str) const;
333 
334 protected:
335 
336 #ifdef WARPX_QED
340  void doQedBreitWheeler (int lev,
341  const amrex::MultiFab& Ex,
342  const amrex::MultiFab& Ey,
343  const amrex::MultiFab& Ez,
344  const amrex::MultiFab& Bx,
345  const amrex::MultiFab& By,
346  const amrex::MultiFab& Bz);
347 
351  void doQedQuantumSync (int lev,
352  const amrex::MultiFab& Ex,
353  const amrex::MultiFab& Ey,
354  const amrex::MultiFab& Ez,
355  const amrex::MultiFab& Bx,
356  const amrex::MultiFab& By,
357  const amrex::MultiFab& Bz);
358 #endif
359 
360  // Particle container types
362 
363  std::vector<std::string> species_names;
364 
365  std::vector<std::string> lasers_names;
366 
367  std::unique_ptr<CollisionHandler> collisionhandler;
368 
370  std::vector<bool> m_deposit_on_main_grid;
371  std::vector<bool> m_laser_deposit_on_main_grid;
372 
374  std::vector<bool> m_gather_from_main_grid;
375 
376  std::vector<PCTypes> species_types;
377 
378  template<typename ...Args>
380  Args const&... pc_dsts) const noexcept
381  {
382  amrex::MFItInfo info;
383 
384  MFItInfoCheckTiling(pc_src, pc_dsts...);
385 
388  }
389 
390 #ifdef AMREX_USE_OMP
391  info.SetDynamic(true);
392 #endif
393 
394  return info;
395  }
396 
397 
398 #ifdef WARPX_QED
399  // The QED engines
400  std::shared_ptr<BreitWheelerEngine> m_shr_p_bw_engine;
401  std::shared_ptr<QuantumSynchrotronEngine> m_shr_p_qs_engine;
402  //_______________________________
403 
408  void InitQED ();
409 
410  //Variables to store how many species need a QED process
413  //________
414 
416  static_cast<amrex::ParticleReal>(
427 
432 
436  void InitQuantumSync ();
437 
441  void InitBreitWheeler ();
442 
448 
454 
475  amrex::Real m_qed_schwinger_xmin = std::numeric_limits<amrex::Real>::lowest();
476  amrex::Real m_qed_schwinger_xmax = std::numeric_limits<amrex::Real>::max();
477  amrex::Real m_qed_schwinger_ymin = std::numeric_limits<amrex::Real>::lowest();
478  amrex::Real m_qed_schwinger_ymax = std::numeric_limits<amrex::Real>::max();
479  amrex::Real m_qed_schwinger_zmin = std::numeric_limits<amrex::Real>::lowest();
480  amrex::Real m_qed_schwinger_zmax = std::numeric_limits<amrex::Real>::max();
481 
482 #endif
483 
484 private:
485 
486  // physical particles (+ laser)
488  // Temporary particle container, used e.g. for particle splitting.
489  std::unique_ptr<PhysicalParticleContainer> pc_tmp;
490 
491  void ReadParameters ();
492 
493  void mapSpeciesProduct ();
494 
496 
497  void MFItInfoCheckTiling(const WarpXParticleContainer& /*pc_src*/) const noexcept
498  {}
499 
500  template<typename First, typename ...Args>
502  First const& pc_dst, Args const&... others) const noexcept
503  {
505  WARPX_ALWAYS_ASSERT_WITH_MESSAGE(pc_dst.do_tiling,
506  "For particle creation processes, either all or none of the "
507  "particle species must use tiling.");
508  }
509 
510  MFItInfoCheckTiling(pc_src, others...);
511  }
512 
519 
520 #ifdef WARPX_QED
526  void CheckQEDProductSpecies();
527 #endif
528 
529 
530 };
531 #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
void RedistributeLocal(int num_ghost)
Definition: MultiParticleContainer.cpp:655
std::unique_ptr< amrex::Parser > m_Bz_particle_parser
Definition: MultiParticleContainer.H:303
void InitQED()
Definition: MultiParticleContainer.cpp:979
std::unique_ptr< amrex::Parser > m_Ey_particle_parser
Definition: MultiParticleContainer.H:306
std::unique_ptr< amrex::MultiFab > GetChargeDensity(int lev, bool local=false)
Definition: MultiParticleContainer.cpp:605
int nLasers() const
Definition: MultiParticleContainer.H:238
amrex::Real m_qed_schwinger_zmax
Definition: MultiParticleContainer.H:480
amrex::Real m_qed_schwinger_zmin
Definition: MultiParticleContainer.H:479
void InitMultiPhysicsModules()
Definition: MultiParticleContainer.cpp:453
std::string m_B_ext_particle_s
Definition: MultiParticleContainer.H:295
void QuantumSyncGenerateTable()
Definition: MultiParticleContainer.cpp:1139
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_strengths_E
Definition: MultiParticleContainer.H:312
void AllocData()
Definition: MultiParticleContainer.cpp:421
void SetParticleDistributionMap(int lev, amrex::DistributionMapping &new_dm)
Definition: MultiParticleContainer.cpp:719
void InitData()
Definition: MultiParticleContainer.cpp:430
amrex::Real m_qed_schwinger_ymin
Definition: MultiParticleContainer.H:477
std::vector< bool > m_laser_deposit_on_main_grid
Definition: MultiParticleContainer.H:371
void Increment(amrex::MultiFab &mf, int lev)
Definition: MultiParticleContainer.cpp:703
void mapSpeciesProduct()
Definition: MultiParticleContainer.cpp:780
void ContinuousInjection(const amrex::RealBox &injection_box) const
Definition: MultiParticleContainer.cpp:732
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_lengths
Definition: MultiParticleContainer.H:315
MultiParticleContainer(MultiParticleContainer const &)=delete
int NSpeciesBreitWheeler() const
Definition: MultiParticleContainer.H:431
void ReadHeader(std::istream &is)
Definition: ParticleIO.cpp:218
static constexpr auto m_default_quantum_sync_photon_creation_energy_threshold
Definition: MultiParticleContainer.H:415
void ApplyBoundaryConditions()
Definition: MultiParticleContainer.cpp:663
std::vector< std::string > lasers_names
Definition: MultiParticleContainer.H:365
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_starts
Definition: MultiParticleContainer.H:314
bool m_do_qed_schwinger
Definition: MultiParticleContainer.H:456
std::vector< std::string > GetLasersNames() const
Definition: MultiParticleContainer.H:283
std::unique_ptr< amrex::Parser > m_By_particle_parser
Definition: MultiParticleContainer.H:302
void PushX(amrex::Real dt)
Definition: MultiParticleContainer.cpp:498
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_lengths
Definition: MultiParticleContainer.H:311
std::string m_qed_schwinger_ele_product_name
Definition: MultiParticleContainer.H:458
int m_nspecies_quantum_sync
Definition: MultiParticleContainer.H:411
int m_qed_schwinger_ele_product
Definition: MultiParticleContainer.H:462
int nSpecies() const
Definition: MultiParticleContainer.H:237
int m_nspecies_breit_wheeler
Definition: MultiParticleContainer.H:412
void MFItInfoCheckTiling(const WarpXParticleContainer &) const noexcept
Definition: MultiParticleContainer.H:497
void InitQuantumSync()
Definition: MultiParticleContainer.cpp:1008
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:540
std::shared_ptr< BreitWheelerEngine > m_shr_p_bw_engine
Definition: MultiParticleContainer.H:400
amrex::Vector< amrex::ParticleReal > m_E_external_particle
Definition: MultiParticleContainer.H:299
void Redistribute()
Definition: MultiParticleContainer.cpp:639
std::vector< PCTypes > species_types
Definition: MultiParticleContainer.H:376
void PostRestart()
Definition: MultiParticleContainer.cpp:442
int NSpeciesQuantumSync() const
Definition: MultiParticleContainer.H:426
amrex::ParticleReal m_quantum_sync_photon_creation_energy_threshold
Definition: MultiParticleContainer.H:420
WarpXParticleContainer & GetParticleContainerFromName(const std::string &name) const
Definition: MultiParticleContainer.cpp:410
MultiParticleContainer & operator=(MultiParticleContainer &&)=default
std::vector< std::string > species_names
Definition: MultiParticleContainer.H:363
amrex::Box ComputeSchwingerGlobalBox() const
Definition: MultiParticleContainer.cpp:1419
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_strengths_B
Definition: MultiParticleContainer.H:317
amrex::Real m_qed_schwinger_ymax
Definition: MultiParticleContainer.H:478
void doQEDSchwinger()
Definition: MultiParticleContainer.cpp:1314
std::string m_E_ext_particle_s
Definition: MultiParticleContainer.H:296
amrex::Real m_qed_schwinger_xmin
Definition: MultiParticleContainer.H:475
int nSpeciesGatherFromMainGrid() const
Definition: MultiParticleContainer.H:259
void defineAllParticleTiles()
Definition: MultiParticleContainer.cpp:647
void BreitWheelerGenerateTable()
Definition: MultiParticleContainer.cpp:1229
std::unique_ptr< PhysicalParticleContainer > pc_tmp
Definition: MultiParticleContainer.H:489
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:516
void ReadParameters()
Definition: MultiParticleContainer.cpp:129
int doContinuousInjection() const
Definition: MultiParticleContainer.cpp:752
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_strengths_B
Definition: MultiParticleContainer.H:313
int nSpeciesDepositOnMainGrid() const
Definition: MultiParticleContainer.H:253
void DepositCharge(amrex::Vector< std::unique_ptr< amrex::MultiFab > > &rho, amrex::Real relative_time)
Deposit charge density.
Definition: MultiParticleContainer.cpp:568
std::unique_ptr< CollisionHandler > collisionhandler
Definition: MultiParticleContainer.H:367
amrex::ParticleReal m_repeated_plasma_lens_period
Definition: MultiParticleContainer.H:309
amrex::Real m_qed_schwinger_xmax
Definition: MultiParticleContainer.H:476
WarpXParticleContainer * GetParticleContainerPtr(int index) const
Definition: MultiParticleContainer.H:82
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:506
amrex::Vector< amrex::ParticleReal > m_B_external_particle
Definition: MultiParticleContainer.H:298
void MFItInfoCheckTiling(const WarpXParticleContainer &pc_src, First const &pc_dst, Args const &... others) const noexcept
Definition: MultiParticleContainer.H:501
int getSpeciesID(std::string product_str) const
Definition: MultiParticleContainer.cpp:825
std::vector< std::string > GetSpeciesAndLasersNames() const
Definition: MultiParticleContainer.H:285
~MultiParticleContainer()
Definition: MultiParticleContainer.H:71
amrex::Vector< amrex::Long > NumberOfParticlesInGrid(int lev) const
Definition: MultiParticleContainer.cpp:680
std::vector< std::string > GetSpeciesNames() const
Definition: MultiParticleContainer.H:281
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:873
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:671
bool m_do_back_transformed_particles
Definition: MultiParticleContainer.H:495
void SortParticlesByBin(amrex::IntVect bin_size)
Definition: MultiParticleContainer.cpp:627
amrex::Vector< std::unique_ptr< WarpXParticleContainer > > allcontainers
Definition: MultiParticleContainer.H:487
void SetParticleBoxArray(int lev, amrex::BoxArray &new_ba)
Definition: MultiParticleContainer.cpp:711
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:1576
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:742
MultiParticleContainer(MultiParticleContainer &&)=default
void doResampling(int timestep, bool verbose)
This function loops over all species and performs resampling if appropriate.
Definition: MultiParticleContainer.cpp:945
void CheckIonizationProductSpecies()
Definition: MultiParticleContainer.cpp:956
std::unique_ptr< amrex::Parser > m_Ex_particle_parser
Definition: MultiParticleContainer.H:305
std::array< amrex::ParticleReal, 3 > meanParticleVelocity(int index)
Definition: MultiParticleContainer.H:87
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_starts
Definition: MultiParticleContainer.H:310
int m_qed_schwinger_threshold_poisson_gaussian
Definition: MultiParticleContainer.H:471
std::vector< bool > m_gather_from_main_grid
instead of gathering fields from the finest patch level, gather from the coarsest
Definition: MultiParticleContainer.H:374
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:1481
void CheckQEDProductSpecies()
Definition: MultiParticleContainer.cpp:1655
MultiParticleContainer & operator=(MultiParticleContainer const &)=delete
amrex::Real m_qed_schwinger_y_size
Definition: MultiParticleContainer.H:466
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_strengths_E
Definition: MultiParticleContainer.H:316
void ScrapeParticles(const amrex::Vector< const amrex::MultiFab * > &distance_to_eb)
Definition: MultiParticleContainer.cpp:967
PhysicalParticleContainer & GetPCtmp()
Definition: MultiParticleContainer.H:291
int nContainers() const
Definition: MultiParticleContainer.H:239
void WriteHeader(std::ostream &os) const
Definition: ParticleIO.cpp:229
std::string m_qed_schwinger_pos_product_name
Definition: MultiParticleContainer.H:460
void Restart(const std::string &dir)
Definition: ParticleIO.cpp:122
void ContinuousFluxInjection(amrex::Real t, amrex::Real dt) const
Definition: MultiParticleContainer.cpp:768
int m_qed_schwinger_pos_product
Definition: MultiParticleContainer.H:464
amrex::MFItInfo getMFItInfo(const WarpXParticleContainer &pc_src, Args const &... pc_dsts) const noexcept
Definition: MultiParticleContainer.H:379
MultiParticleContainer(amrex::AmrCore *amr_core)
Definition: MultiParticleContainer.cpp:93
WarpXParticleContainer & GetParticleContainer(int index) const
Definition: MultiParticleContainer.H:79
void doCollisions(amrex::Real cur_time, amrex::Real dt)
Definition: MultiParticleContainer.cpp:939
std::unique_ptr< amrex::Parser > m_Ez_particle_parser
Definition: MultiParticleContainer.H:307
std::shared_ptr< QuantumSynchrotronEngine > m_shr_p_qs_engine
Definition: MultiParticleContainer.H:401
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:370
void SetDoBackTransformedParticles(bool do_back_transformed_particles)
Definition: MultiParticleContainer.cpp:849
PCTypes
Definition: MultiParticleContainer.H:361
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:471
std::unique_ptr< amrex::Parser > m_Bx_particle_parser
Definition: MultiParticleContainer.H:301
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:1495
Definition: PhysicalParticleContainer.H:46
Definition: WarpXParticleContainer.H:104
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