WarpX
Loading...
Searching...
No Matches
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
16
18#ifdef WARPX_QED
21#endif
23#include "Utils/TextMsg.H"
24#include "Utils/WarpXConst.H"
26#include "ParticleBoundaries.H"
28
30
31#include <AMReX_BLassert.H>
32#include <AMReX_Box.H>
33#include <AMReX_Config.H>
34#include <AMReX_Geometry.H>
35#include <AMReX_GpuControl.H>
36#include <AMReX_INT.H>
37#include <AMReX_MFIter.H>
38#include <AMReX_REAL.H>
39#include <AMReX_RealBox.H>
40#include <AMReX_Vector.H>
41
42#include <AMReX_BaseFwd.H>
43#include <AMReX_AmrCoreFwd.H>
44
45#include <algorithm>
46#include <array>
47#include <iosfwd>
48#include <iterator>
49#include <limits>
50#include <memory>
51#include <string>
52#include <vector>
53
69{
70
71public:
72
74
76
81
82 [[nodiscard]] WarpXParticleContainer&
83 GetParticleContainer (int index) const {return *allcontainers[index];}
84
85 [[nodiscard]] WarpXParticleContainer*
86 GetParticleContainerPtr (int index) const {return allcontainers[index].get();}
87
88 [[nodiscard]] WarpXParticleContainer&
89 GetParticleContainerFromName (const std::string& name) const;
90
91 std::array<amrex::ParticleReal, 3> meanParticleVelocity(int index) {
92 return allcontainers[index]->meanParticleVelocity();
93 }
94
96
97 void AllocData ();
98
99 void InitData ();
100
102
108 void Evolve (
110 int lev,
111 std::string const& current_fp_string,
112 amrex::Real t,
113 amrex::Real dt,
114 SubcyclingHalf subcycling_half=SubcyclingHalf::None,
115 bool skip_deposition=false,
116 PositionPushType position_push_type=PositionPushType::Full,
117 MomentumPushType momentum_push_type=MomentumPushType::Full,
118 ImplicitOptions const * implicit_options = nullptr
119 );
120
125 void PushX (amrex::Real dt);
126
133 void PushP (int lev, amrex::Real dt,
134 const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
135 const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz);
136
143 std::unique_ptr<amrex::MultiFab> GetZeroChargeDensity(int lev);
144
154 void
156 amrex::Real relative_time);
157
169 void
171 amrex::Real dt, amrex::Real relative_time);
172
183 void
185
189
190 std::unique_ptr<amrex::MultiFab> GetChargeDensity(int lev, bool local = false);
191
193
194 void doFieldIonization (int lev,
195 const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
196 const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz);
197
198 void doCollisions (int step, amrex::Real cur_time, amrex::Real dt);
199
207 void doResampling (const amrex::Vector<amrex::Geometry>& geom, int timestep, bool verbose);
208
209#ifdef WARPX_QED
217 void doQEDSchwinger ();
218
223 [[nodiscard]] amrex::Box ComputeSchwingerGlobalBox () const;
224#endif
225
226 void Restart (const std::string& dir);
227
228 void PostRestart ();
229
230 void ReadHeader (std::istream& is);
231
232 void WriteHeader (std::ostream& os) const;
233
234 void SortParticlesByBin (
235 const amrex::IntVect& bin_size,
236 bool sort_particles_for_deposition,
237 const amrex::IntVect& sort_idx_type);
238
239 void Redistribute ();
240
242
244
245 void RedistributeLocal (int num_ghost);
246
250
258 [[nodiscard]] amrex::Vector<amrex::Long> GetZeroParticlesInGrid(int lev) const;
259
260 [[nodiscard]] amrex::Vector<amrex::Long> NumberOfParticlesInGrid(int lev) const;
261
262 void Increment (amrex::MultiFab& mf, int lev);
263
264 void SetParticleBoxArray (int lev, amrex::BoxArray& new_ba);
266
267 [[nodiscard]] int nSpecies () const {return static_cast<int>(species_names.size());}
268 [[nodiscard]] int nLasers () const {return static_cast<int>(lasers_names.size());}
269 [[nodiscard]] int nContainers () const {return static_cast<int>(allcontainers.size());}
270
275 void SetDoBackTransformedParticles (bool do_back_transformed_particles);
281 void SetDoBackTransformedParticles (const std::string& species_name, bool do_back_transformed_particles);
282
285 [[nodiscard]] bool getDoBackTransformedParticles () const
286 {
288 }
289
290 [[nodiscard]] int nSpeciesDepositOnMainGrid () const
291 {
292 bool const onMainGrid = true;
293 auto const & v = m_deposit_on_main_grid;
294 return static_cast<int>(std::count( v.begin(), v.end(), onMainGrid ));
295 }
296
297 [[nodiscard]] int nSpeciesGatherFromMainGrid() const
298 {
299 bool const fromMainGrid = true;
300 auto const & v = m_gather_from_main_grid;
301 return static_cast<int>(std::count( v.begin(), v.end(), fromMainGrid ));
302 }
303
304 // Inject particles during the simulation (for particles entering the
305 // simulation domain after some iterations, due to flowing plasma and/or
306 // moving window).
307 void ContinuousInjection(const amrex::RealBox& injection_box) const;
308
313 void UpdateAntennaPosition(amrex::Real dt) const;
314
315 [[nodiscard]] int doContinuousInjection() const;
316
317 // Inject particles from a surface during the simulation
319
320 [[nodiscard]] std::vector<std::string> GetSpeciesNames() const { return species_names; }
321
322 [[nodiscard]] std::vector<std::string> GetLasersNames() const { return lasers_names; }
323
324 [[nodiscard]] std::vector<std::string> GetSpeciesAndLasersNames() const
325 {
326 std::vector<std::string> tmp = species_names;
327 tmp.insert(tmp.end(), lasers_names.begin(), lasers_names.end());
328 return tmp;
329 }
330
332
333 std::string m_B_ext_particle_s = "none";
334 std::string m_E_ext_particle_s = "none";
335 // Parser for B_external on the particle
336 std::unique_ptr<amrex::Parser> m_Bx_particle_parser;
337 std::unique_ptr<amrex::Parser> m_By_particle_parser;
338 std::unique_ptr<amrex::Parser> m_Bz_particle_parser;
339 // Parser for E_external on the particle
340 std::unique_ptr<amrex::Parser> m_Ex_particle_parser;
341 std::unique_ptr<amrex::Parser> m_Ey_particle_parser;
342 std::unique_ptr<amrex::Parser> m_Ez_particle_parser;
343
345
355
356#ifdef WARPX_QED
360 void doQedEvents (int lev,
361 const amrex::MultiFab& Ex,
362 const amrex::MultiFab& Ey,
363 const amrex::MultiFab& Ez,
364 const amrex::MultiFab& Bx,
365 const amrex::MultiFab& By,
366 const amrex::MultiFab& Bz);
367#endif
368
369 [[nodiscard]] int getSpeciesID (const std::string& product_str) const;
370
373
374protected:
375
376#ifdef WARPX_QED
380 void doQedBreitWheeler (int lev,
381 const amrex::MultiFab& Ex,
382 const amrex::MultiFab& Ey,
383 const amrex::MultiFab& Ez,
384 const amrex::MultiFab& Bx,
385 const amrex::MultiFab& By,
386 const amrex::MultiFab& Bz);
387
391 void doQedQuantumSync (int lev,
392 const amrex::MultiFab& Ex,
393 const amrex::MultiFab& Ey,
394 const amrex::MultiFab& Ez,
395 const amrex::MultiFab& Bx,
396 const amrex::MultiFab& By,
397 const amrex::MultiFab& Bz);
398#endif
399
400 // Particle container types
402
403 std::vector<std::string> species_names;
404
405 std::vector<std::string> lasers_names;
406
407 std::unique_ptr<CollisionHandler> collisionhandler;
408
410 std::vector<bool> m_deposit_on_main_grid;
412
414 std::vector<bool> m_gather_from_main_grid;
415
416 std::vector<PCTypes> species_types;
417
418 template<typename ...Args>
419 [[nodiscard]]
421 Args const&... pc_dsts) const noexcept
422 {
423 amrex::MFItInfo info;
424
425 MFItInfoCheckTiling(pc_src, pc_dsts...);
426
429 }
430
431#ifdef AMREX_USE_OMP
432 info.SetDynamic(true);
433#endif
434
435 return info;
436 }
437
438
439#ifdef WARPX_QED
440 // The QED engines
441 std::shared_ptr<BreitWheelerEngine> m_shr_p_bw_engine;
442 std::shared_ptr<QuantumSynchrotronEngine> m_shr_p_qs_engine;
443 //_______________________________
444
449 void InitQED ();
450
451 //Variables to store how many species need a QED process
454 //________
455
457 static_cast<amrex::ParticleReal>(
459
460
463
467 [[nodiscard]] int NSpeciesQuantumSync() const { return m_nspecies_quantum_sync;}
468
472 [[nodiscard]] int NSpeciesBreitWheeler() const { return m_nspecies_breit_wheeler;}
473
477 void InitQuantumSync ();
478
482 void InitBreitWheeler ();
483
489
495
497 bool m_do_qed_schwinger = false;
516 amrex::Real m_qed_schwinger_xmin = std::numeric_limits<amrex::Real>::lowest();
517 amrex::Real m_qed_schwinger_xmax = std::numeric_limits<amrex::Real>::max();
518 amrex::Real m_qed_schwinger_ymin = std::numeric_limits<amrex::Real>::lowest();
519 amrex::Real m_qed_schwinger_ymax = std::numeric_limits<amrex::Real>::max();
520 amrex::Real m_qed_schwinger_zmin = std::numeric_limits<amrex::Real>::lowest();
521 amrex::Real m_qed_schwinger_zmax = std::numeric_limits<amrex::Real>::max();
522
523#endif
524
525private:
526
527 // physical particles (+ laser)
529
530 void ReadParameters ();
531
532 void mapSpeciesProduct ();
533
535
536 void MFItInfoCheckTiling(const WarpXParticleContainer& /*pc_src*/) const noexcept
537 {}
538
539 template<typename First, typename ...Args>
541 First const& pc_dst, Args const&... others) const noexcept
542 {
544 WARPX_ALWAYS_ASSERT_WITH_MESSAGE(pc_dst.do_tiling,
545 "For particle creation processes, either all or none of the "
546 "particle species must use tiling.");
547 }
548
549 MFItInfoCheckTiling(pc_src, others...);
550 }
551
558
559#ifdef WARPX_QED
566#endif
567
568
569};
570#endif /*WARPX_ParticleContainer_H_*/
@ v
Definition RigidInjectedParticleContainer.H:27
#define WARPX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition TextMsg.H:13
PositionPushType
For advanced collision algorithms that split the particle push in substeps.
Definition WarpXAlgorithmSelection.H:181
@ Full
Definition WarpXAlgorithmSelection.H:181
MomentumPushType
For advanced collision algorithms that split the particle push in substeps.
Definition WarpXAlgorithmSelection.H:189
@ Full
Definition WarpXAlgorithmSelection.H:189
SubcyclingHalf
Subcycling half selector.
Definition WarpXAlgorithmSelection.H:166
@ None
Definition WarpXAlgorithmSelection.H:166
This class contains the parameters for the external particle fields.
Definition ExternalParticleFields.H:40
void SortParticlesByBin(const amrex::IntVect &bin_size, bool sort_particles_for_deposition, const amrex::IntVect &sort_idx_type)
Definition MultiParticleContainer.cpp:750
std::vector< std::string > GetLasersNames() const
Definition MultiParticleContainer.H:322
void RedistributeLocal(int num_ghost)
Definition MultiParticleContainer.cpp:790
amrex::Vector< std::unique_ptr< WarpXParticleContainer > >::iterator begin()
Definition MultiParticleContainer.H:371
void deleteInvalidParticles()
Definition MultiParticleContainer.cpp:782
std::unique_ptr< amrex::Parser > m_Bz_particle_parser
Definition MultiParticleContainer.H:338
void InitQED()
Definition MultiParticleContainer.cpp:1135
std::unique_ptr< amrex::Parser > m_Ey_particle_parser
Definition MultiParticleContainer.H:341
std::unique_ptr< amrex::MultiFab > GetChargeDensity(int lev, bool local=false)
Definition MultiParticleContainer.cpp:653
int nLasers() const
Definition MultiParticleContainer.H:268
bool getDoBackTransformedParticles() const
Definition MultiParticleContainer.H:285
amrex::Real m_qed_schwinger_zmax
Definition MultiParticleContainer.H:521
amrex::Real m_qed_schwinger_zmin
Definition MultiParticleContainer.H:520
void InitMultiPhysicsModules()
Definition MultiParticleContainer.cpp:452
std::string m_B_ext_particle_s
Definition MultiParticleContainer.H:333
void QuantumSyncGenerateTable()
Definition MultiParticleContainer.cpp:1298
MultiParticleContainer & operator=(MultiParticleContainer const &)=delete
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_strengths_E
Definition MultiParticleContainer.H:349
void AllocData()
Definition MultiParticleContainer.cpp:424
void SetParticleDistributionMap(int lev, amrex::DistributionMapping &new_dm)
Definition MultiParticleContainer.cpp:854
void InitData()
Definition MultiParticleContainer.cpp:432
amrex::Real m_qed_schwinger_ymin
Definition MultiParticleContainer.H:518
std::vector< bool > m_laser_deposit_on_main_grid
Definition MultiParticleContainer.H:411
void Evolve(ablastr::fields::MultiFabRegister &fields, int lev, std::string const &current_fp_string, amrex::Real t, amrex::Real dt, SubcyclingHalf subcycling_half=SubcyclingHalf::None, bool skip_deposition=false, PositionPushType position_push_type=PositionPushType::Full, MomentumPushType momentum_push_type=MomentumPushType::Full, ImplicitOptions const *implicit_options=nullptr)
This evolves all the particles by one PIC time step, including current deposition,...
Definition MultiParticleContainer.cpp:470
void Increment(amrex::MultiFab &mf, int lev)
Definition MultiParticleContainer.cpp:838
void mapSpeciesProduct()
Definition MultiParticleContainer.cpp:915
void ContinuousInjection(const amrex::RealBox &injection_box) const
Definition MultiParticleContainer.cpp:867
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_lengths
Definition MultiParticleContainer.H:352
MultiParticleContainer(MultiParticleContainer const &)=delete
void doResampling(const amrex::Vector< amrex::Geometry > &geom, int timestep, bool verbose)
This function loops over all species and performs resampling if appropriate.
Definition MultiParticleContainer.cpp:1103
amrex::Vector< std::unique_ptr< WarpXParticleContainer > >::iterator end()
Definition MultiParticleContainer.H:372
int NSpeciesBreitWheeler() const
Definition MultiParticleContainer.H:472
void ReadHeader(std::istream &is)
Definition ParticleIO.cpp:231
static constexpr auto m_default_quantum_sync_photon_creation_energy_threshold
Definition MultiParticleContainer.H:456
void ApplyBoundaryConditions()
Definition MultiParticleContainer.cpp:798
std::vector< std::string > lasers_names
Definition MultiParticleContainer.H:405
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_starts
Definition MultiParticleContainer.H:351
bool m_do_qed_schwinger
Definition MultiParticleContainer.H:497
std::unique_ptr< amrex::Parser > m_By_particle_parser
Definition MultiParticleContainer.H:337
void PushX(amrex::Real dt)
This pushes the particle positions by one time step for all the species in the MultiParticleContainer...
Definition MultiParticleContainer.cpp:515
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_lengths
Definition MultiParticleContainer.H:348
std::string m_qed_schwinger_ele_product_name
Definition MultiParticleContainer.H:499
int m_nspecies_quantum_sync
Definition MultiParticleContainer.H:452
int m_qed_schwinger_ele_product
Definition MultiParticleContainer.H:503
int nSpecies() const
Definition MultiParticleContainer.H:267
int m_nspecies_breit_wheeler
Definition MultiParticleContainer.H:453
void MFItInfoCheckTiling(const WarpXParticleContainer &) const noexcept
Definition MultiParticleContainer.H:536
void InitQuantumSync()
Definition MultiParticleContainer.cpp:1166
std::shared_ptr< BreitWheelerEngine > m_shr_p_bw_engine
Definition MultiParticleContainer.H:441
void Redistribute()
Definition MultiParticleContainer.cpp:766
std::vector< PCTypes > species_types
Definition MultiParticleContainer.H:416
void PostRestart()
Definition MultiParticleContainer.cpp:442
int NSpeciesQuantumSync() const
Definition MultiParticleContainer.H:467
amrex::ParticleReal m_quantum_sync_photon_creation_energy_threshold
Definition MultiParticleContainer.H:461
void DepositCurrent(ablastr::fields::MultiLevelVectorField const &J, amrex::Real dt, amrex::Real relative_time)
Deposit current density.
Definition MultiParticleContainer.cpp:558
WarpXParticleContainer & GetParticleContainerFromName(const std::string &name) const
Definition MultiParticleContainer.cpp:404
std::vector< std::string > species_names
Definition MultiParticleContainer.H:403
amrex::Box ComputeSchwingerGlobalBox() const
Definition MultiParticleContainer.cpp:1579
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_strengths_B
Definition MultiParticleContainer.H:354
amrex::Real m_qed_schwinger_ymax
Definition MultiParticleContainer.H:519
void doQEDSchwinger()
Definition MultiParticleContainer.cpp:1471
std::string m_E_ext_particle_s
Definition MultiParticleContainer.H:334
amrex::Real m_qed_schwinger_xmin
Definition MultiParticleContainer.H:516
int nSpeciesGatherFromMainGrid() const
Definition MultiParticleContainer.H:297
WarpXParticleContainer * GetParticleContainerPtr(int index) const
Definition MultiParticleContainer.H:86
void defineAllParticleTiles()
Definition MultiParticleContainer.cpp:774
void BreitWheelerGenerateTable()
Definition MultiParticleContainer.cpp:1387
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:533
void ReadParameters()
Definition MultiParticleContainer.cpp:130
int doContinuousInjection() const
Definition MultiParticleContainer.cpp:887
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_strengths_B
Definition MultiParticleContainer.H:350
int nSpeciesDepositOnMainGrid() const
Definition MultiParticleContainer.H:290
void DepositTemperatures(ablastr::fields::MultiFabRegister &fields, amrex::Real relative_time)
Deposit temperature to species MFs. This is done for each species and can be used in the future for a...
Definition MultiParticleContainer.cpp:623
MultiParticleContainer & operator=(MultiParticleContainer &&)=default
std::unique_ptr< CollisionHandler > collisionhandler
Definition MultiParticleContainer.H:407
amrex::ParticleReal m_repeated_plasma_lens_period
Definition MultiParticleContainer.H:346
amrex::Real m_qed_schwinger_xmax
Definition MultiParticleContainer.H:517
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:523
void doCollisions(int step, amrex::Real cur_time, amrex::Real dt)
Definition MultiParticleContainer.cpp:1097
void MFItInfoCheckTiling(const WarpXParticleContainer &pc_src, First const &pc_dst, Args const &... others) const noexcept
Definition MultiParticleContainer.H:540
int getSpeciesID(const std::string &product_str) const
Definition MultiParticleContainer.cpp:966
amrex::Vector< amrex::Long > NumberOfParticlesInGrid(int lev) const
Definition MultiParticleContainer.cpp:815
void InitBreitWheeler()
Definition MultiParticleContainer.cpp:1238
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:1031
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:806
bool m_do_back_transformed_particles
Definition MultiParticleContainer.H:534
void GenerateGlobalDebyeLength()
Definition MultiParticleContainer.cpp:675
std::array< amrex::ParticleReal, 3 > meanParticleVelocity(int index)
Definition MultiParticleContainer.H:91
amrex::Vector< std::unique_ptr< WarpXParticleContainer > > allcontainers
Definition MultiParticleContainer.H:528
void SetParticleBoxArray(int lev, amrex::BoxArray &new_ba)
Definition MultiParticleContainer.cpp:846
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:1738
WarpXParticleContainer & GetParticleContainer(int index) const
Definition MultiParticleContainer.H:83
void ScrapeParticlesAtEB(ablastr::fields::MultiLevelScalarField const &distance_to_eb)
Definition MultiParticleContainer.cpp:1126
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:877
MultiParticleContainer(MultiParticleContainer &&)=default
void CheckIonizationProductSpecies()
Definition MultiParticleContainer.cpp:1115
std::unique_ptr< amrex::Parser > m_Ex_particle_parser
Definition MultiParticleContainer.H:340
std::vector< std::string > GetSpeciesAndLasersNames() const
Definition MultiParticleContainer.H:324
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_starts
Definition MultiParticleContainer.H:347
int m_qed_schwinger_threshold_poisson_gaussian
Definition MultiParticleContainer.H:512
std::vector< bool > m_gather_from_main_grid
instead of gathering fields from the finest patch level, gather from the coarsest
Definition MultiParticleContainer.H:414
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:1641
void DepositCharge(const ablastr::fields::MultiLevelScalarField &rho, amrex::Real relative_time)
Deposit charge density.
Definition MultiParticleContainer.cpp:586
void CheckQEDProductSpecies()
Definition MultiParticleContainer.cpp:1819
amrex::Real m_qed_schwinger_y_size
Definition MultiParticleContainer.H:507
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_strengths_E
Definition MultiParticleContainer.H:353
amrex::ParticleReal maxParticleVelocity()
Definition MultiParticleContainer.cpp:415
int nContainers() const
Definition MultiParticleContainer.H:269
ExternalParticleFields m_external_particle_fields_metadata
Definition MultiParticleContainer.H:344
void WriteHeader(std::ostream &os) const
Definition ParticleIO.cpp:242
std::string m_qed_schwinger_pos_product_name
Definition MultiParticleContainer.H:501
void Restart(const std::string &dir)
Definition ParticleIO.cpp:129
std::vector< std::string > GetSpeciesNames() const
Definition MultiParticleContainer.H:320
void ContinuousFluxInjection(amrex::Real t, amrex::Real dt) const
Definition MultiParticleContainer.cpp:903
int m_qed_schwinger_pos_product
Definition MultiParticleContainer.H:505
amrex::MFItInfo getMFItInfo(const WarpXParticleContainer &pc_src, Args const &... pc_dsts) const noexcept
Definition MultiParticleContainer.H:420
MultiParticleContainer(amrex::AmrCore *amr_core)
Definition MultiParticleContainer.cpp:96
std::unique_ptr< amrex::Parser > m_Ez_particle_parser
Definition MultiParticleContainer.H:342
std::shared_ptr< QuantumSynchrotronEngine > m_shr_p_qs_engine
Definition MultiParticleContainer.H:442
~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:410
void SetDoBackTransformedParticles(bool do_back_transformed_particles)
Definition MultiParticleContainer.cpp:990
PCTypes
Definition MultiParticleContainer.H:401
@ RigidInjected
Definition MultiParticleContainer.H:401
@ Physical
Definition MultiParticleContainer.H:401
@ Photon
Definition MultiParticleContainer.H:401
std::unique_ptr< amrex::Parser > m_Bx_particle_parser
Definition MultiParticleContainer.H:336
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:1655
Definition WarpXParticleContainer.H:195
amrex_real Real
amrex_particle_real ParticleReal
PODVector< T, ArenaAllocator< T > > DeviceVector
constexpr auto c
vacuum speed of light [m/s]
Definition constant.H:153
constexpr auto m_e
electron mass [kg]
Definition constant.H:165
amrex::Vector< ScalarField > MultiLevelScalarField
Definition MultiFabRegister.H:200
amrex::Vector< VectorField > MultiLevelVectorField
Definition MultiFabRegister.H:208
bool notInLaunchRegion() noexcept
BoxND< 3 > Box
IntVectND< 3 > IntVect
Definition ImplicitOptions.H:7
Definition MultiFabRegister.H:262
MFItInfo & SetDynamic(bool f) noexcept
MFItInfo & EnableTiling(const IntVect &ts=FabArrayBase::mfiter_tile_size) noexcept