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
73 explicit MultiParticleContainer (amrex::AmrCore* amr_core);
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 TransformMomentumToCurvilinear (bool forward);
98
99 void AllocData ();
100
101 void InitData ();
102
104
110 void Evolve (
112 int lev,
113 std::string const& current_fp_string,
114 amrex::Real t,
115 amrex::Real dt,
116 SubcyclingHalf subcycling_half=SubcyclingHalf::None,
117 bool skip_deposition=false,
118 PositionPushType position_push_type=PositionPushType::Full,
119 MomentumPushType momentum_push_type=MomentumPushType::Full,
120 ImplicitOptions const * implicit_options = nullptr
121 );
122
127 int lev, amrex::Real dt);
128
133 void PushX (amrex::Real dt);
134
141 void PushP (int lev, amrex::Real dt,
142 const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
143 const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz,
144 MomentumPushType momentum_push_type);
145
152 std::unique_ptr<amrex::MultiFab> GetZeroChargeDensity(int lev);
153
163 void
165 amrex::Real relative_time);
166
178 void
180 amrex::Real dt, amrex::Real relative_time);
181
192 void
194
198
199 std::unique_ptr<amrex::MultiFab> GetChargeDensity(int lev, bool local = false);
200
201 std::unique_ptr<amrex::MultiFab> GetGlobalPlasmaFrequency (int lev);
203
204 /*
205 * \brief Calculates the electron-ion collision frequency for the specified species
206 * @param[in] electron_species_name name of the electon species
207 */
208 void CalculateNuei (std::string const & electron_species_name);
209
210 void doFieldIonization (int lev,
211 const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
212 const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz);
213
214 void doCollisions (int step, amrex::Real cur_time, amrex::Real dt);
215
223 void doResampling (const amrex::Vector<amrex::Geometry>& geom, int timestep, bool verbose);
224
225#ifdef WARPX_QED
233 void doQEDSchwinger ();
234
239 [[nodiscard]] amrex::Box ComputeSchwingerGlobalBox () const;
240#endif
241
242 void Restart (const std::string& dir);
243
244 void PostRestart ();
245
246 void ReadHeader (std::istream& is);
247
248 void WriteHeader (std::ostream& os) const;
249
250 void SortParticlesByBin (
251 const amrex::IntVect& bin_size,
252 bool sort_particles_for_deposition,
253 const amrex::IntVect& sort_idx_type);
254
255 void Redistribute ();
256
258
260
261 void RedistributeLocal (int max_cells_travelled);
262
266
274 [[nodiscard]] amrex::Vector<amrex::Long> GetZeroParticlesInGrid(int lev) const;
275
276 [[nodiscard]] amrex::Vector<amrex::Long> NumberOfParticlesInGrid(int lev) const;
277
278 void Increment (amrex::MultiFab& mf, int lev);
279
280 void SetParticleBoxArray (int lev, amrex::BoxArray& new_ba);
282
283 [[nodiscard]] int nSpecies () const {return static_cast<int>(species_names.size());}
284 [[nodiscard]] int nLasers () const {return static_cast<int>(lasers_names.size());}
285 [[nodiscard]] int nContainers () const {return static_cast<int>(allcontainers.size());}
286
291 void SetDoBackTransformedParticles (bool do_back_transformed_particles);
297 void SetDoBackTransformedParticles (const std::string& species_name, bool do_back_transformed_particles);
298
301 [[nodiscard]] bool getDoBackTransformedParticles () const
302 {
304 }
305
306 [[nodiscard]] int nSpeciesDepositOnMainGrid () const
307 {
308 bool const onMainGrid = true;
309 auto const & v = m_deposit_on_main_grid;
310 return static_cast<int>(std::count( v.begin(), v.end(), onMainGrid ));
311 }
312
313 [[nodiscard]] int nSpeciesGatherFromMainGrid() const
314 {
315 bool const fromMainGrid = true;
316 auto const & v = m_gather_from_main_grid;
317 return static_cast<int>(std::count( v.begin(), v.end(), fromMainGrid ));
318 }
319
320 // Inject particles during the simulation (for particles entering the
321 // simulation domain after some iterations, due to flowing plasma and/or
322 // moving window).
323 void ContinuousInjection(const amrex::RealBox& injection_box) const;
324
329 void UpdateAntennaPosition(amrex::Real dt) const;
330
331 [[nodiscard]] int doContinuousInjection() const;
332
333 // Inject particles from a surface during the simulation
335
336 [[nodiscard]] std::vector<std::string> GetSpeciesNames() const { return species_names; }
337
338 [[nodiscard]] std::vector<std::string> GetLasersNames() const { return lasers_names; }
339
340 [[nodiscard]] std::vector<std::string> GetSpeciesAndLasersNames() const
341 {
342 std::vector<std::string> tmp = species_names;
343 tmp.insert(tmp.end(), lasers_names.begin(), lasers_names.end());
344 return tmp;
345 }
346
348
349 std::string m_B_ext_particle_s = "none";
350 std::string m_E_ext_particle_s = "none";
351 // Parser for B_external on the particle
352 std::unique_ptr<amrex::Parser> m_Bx_particle_parser;
353 std::unique_ptr<amrex::Parser> m_By_particle_parser;
354 std::unique_ptr<amrex::Parser> m_Bz_particle_parser;
355 // Parser for E_external on the particle
356 std::unique_ptr<amrex::Parser> m_Ex_particle_parser;
357 std::unique_ptr<amrex::Parser> m_Ey_particle_parser;
358 std::unique_ptr<amrex::Parser> m_Ez_particle_parser;
359
361
371
372#ifdef WARPX_QED
376 void doQedEvents (int lev,
377 const amrex::MultiFab& Ex,
378 const amrex::MultiFab& Ey,
379 const amrex::MultiFab& Ez,
380 const amrex::MultiFab& Bx,
381 const amrex::MultiFab& By,
382 const amrex::MultiFab& Bz);
383#endif
384
385 [[nodiscard]] int getSpeciesID (const std::string& product_str) const;
386
389
390protected:
391
392#ifdef WARPX_QED
396 void doQedBreitWheeler (int lev,
397 const amrex::MultiFab& Ex,
398 const amrex::MultiFab& Ey,
399 const amrex::MultiFab& Ez,
400 const amrex::MultiFab& Bx,
401 const amrex::MultiFab& By,
402 const amrex::MultiFab& Bz);
403
407 void doQedQuantumSync (int lev,
408 const amrex::MultiFab& Ex,
409 const amrex::MultiFab& Ey,
410 const amrex::MultiFab& Ez,
411 const amrex::MultiFab& Bx,
412 const amrex::MultiFab& By,
413 const amrex::MultiFab& Bz);
414#endif
415
416 // Particle container types
418
419 std::vector<std::string> species_names;
420
421 std::vector<std::string> lasers_names;
422
423 std::unique_ptr<CollisionHandler> collisionhandler;
424
426 std::vector<bool> m_deposit_on_main_grid;
428
430 std::vector<bool> m_gather_from_main_grid;
431
432 std::vector<PCTypes> species_types;
433
434 template<typename ...Args>
435 [[nodiscard]]
437 Args const&... pc_dsts) const noexcept
438 {
439 amrex::MFItInfo info;
440
441 MFItInfoCheckTiling(pc_src, pc_dsts...);
442
445 }
446
447#ifdef AMREX_USE_OMP
448 info.SetDynamic(true);
449#endif
450
451 return info;
452 }
453
454
455#ifdef WARPX_QED
456 // The QED engines
457 std::shared_ptr<BreitWheelerEngine> m_shr_p_bw_engine;
458 std::shared_ptr<QuantumSynchrotronEngine> m_shr_p_qs_engine;
459 //_______________________________
460
465 void InitQED ();
466
467 //Variables to store how many species need a QED process
470 //________
471
473 static_cast<amrex::ParticleReal>(
474 2.0 * PhysConst::m_e * PhysConst::c2 );
475
476
479
483 [[nodiscard]] int NSpeciesQuantumSync() const { return m_nspecies_quantum_sync;}
484
488 [[nodiscard]] int NSpeciesBreitWheeler() const { return m_nspecies_breit_wheeler;}
489
493 void InitQuantumSync ();
494
498 void InitBreitWheeler ();
499
505
511
513 bool m_do_qed_schwinger = false;
532 amrex::Real m_qed_schwinger_xmin = std::numeric_limits<amrex::Real>::lowest();
533 amrex::Real m_qed_schwinger_xmax = std::numeric_limits<amrex::Real>::max();
534 amrex::Real m_qed_schwinger_ymin = std::numeric_limits<amrex::Real>::lowest();
535 amrex::Real m_qed_schwinger_ymax = std::numeric_limits<amrex::Real>::max();
536 amrex::Real m_qed_schwinger_zmin = std::numeric_limits<amrex::Real>::lowest();
537 amrex::Real m_qed_schwinger_zmax = std::numeric_limits<amrex::Real>::max();
538
539#endif
540
541private:
542
543 // physical particles (+ laser)
545
546 void ReadParameters ();
547
548 void mapSpeciesProduct ();
549
551
552 void MFItInfoCheckTiling(const WarpXParticleContainer& /*pc_src*/) const noexcept
553 {}
554
555 template<typename First, typename ...Args>
557 First const& pc_dst, Args const&... others) const noexcept
558 {
560 WARPX_ALWAYS_ASSERT_WITH_MESSAGE(pc_dst.do_tiling,
561 "For particle creation processes, either all or none of the "
562 "particle species must use tiling.");
563 }
564
565 MFItInfoCheckTiling(pc_src, others...);
566 }
567
574
575#ifdef WARPX_QED
582#endif
583
584
585};
586#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:180
@ Full
Definition WarpXAlgorithmSelection.H:180
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:861
std::vector< std::string > GetLasersNames() const
Definition MultiParticleContainer.H:338
amrex::Vector< std::unique_ptr< WarpXParticleContainer > >::iterator begin()
Definition MultiParticleContainer.H:387
void deleteInvalidParticles()
Definition MultiParticleContainer.cpp:893
std::unique_ptr< amrex::Parser > m_Bz_particle_parser
Definition MultiParticleContainer.H:354
void InitQED()
Definition MultiParticleContainer.cpp:1248
std::unique_ptr< amrex::Parser > m_Ey_particle_parser
Definition MultiParticleContainer.H:357
std::unique_ptr< amrex::MultiFab > GetChargeDensity(int lev, bool local=false)
Definition MultiParticleContainer.cpp:682
int nLasers() const
Definition MultiParticleContainer.H:284
bool getDoBackTransformedParticles() const
Definition MultiParticleContainer.H:301
amrex::Real m_qed_schwinger_zmax
Definition MultiParticleContainer.H:537
amrex::Real m_qed_schwinger_zmin
Definition MultiParticleContainer.H:536
void InitMultiPhysicsModules()
Definition MultiParticleContainer.cpp:460
std::string m_B_ext_particle_s
Definition MultiParticleContainer.H:349
void QuantumSyncGenerateTable()
Definition MultiParticleContainer.cpp:1411
MultiParticleContainer & operator=(MultiParticleContainer const &)=delete
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_strengths_E
Definition MultiParticleContainer.H:365
void AllocData()
Definition MultiParticleContainer.cpp:432
void SetParticleDistributionMap(int lev, amrex::DistributionMapping &new_dm)
Definition MultiParticleContainer.cpp:967
void InitData()
Definition MultiParticleContainer.cpp:440
amrex::Real m_qed_schwinger_ymin
Definition MultiParticleContainer.H:534
std::vector< bool > m_laser_deposit_on_main_grid
Definition MultiParticleContainer.H:427
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:478
void Increment(amrex::MultiFab &mf, int lev)
Definition MultiParticleContainer.cpp:951
void mapSpeciesProduct()
Definition MultiParticleContainer.cpp:1028
void ContinuousInjection(const amrex::RealBox &injection_box) const
Definition MultiParticleContainer.cpp:980
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_lengths
Definition MultiParticleContainer.H:368
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:1216
void DepositMassMatrices(ablastr::fields::MultiFabRegister &fields, int lev, amrex::Real dt)
Deposit mass matrices.
Definition MultiParticleContainer.cpp:526
amrex::Vector< std::unique_ptr< WarpXParticleContainer > >::iterator end()
Definition MultiParticleContainer.H:388
int NSpeciesBreitWheeler() const
Definition MultiParticleContainer.H:488
void CalculateNuei(std::string const &electron_species_name)
Definition MultiParticleContainer.cpp:820
void ReadHeader(std::istream &is)
Definition ParticleIO.cpp:231
static constexpr auto m_default_quantum_sync_photon_creation_energy_threshold
Definition MultiParticleContainer.H:472
void ApplyBoundaryConditions()
Definition MultiParticleContainer.cpp:911
std::vector< std::string > lasers_names
Definition MultiParticleContainer.H:421
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_starts
Definition MultiParticleContainer.H:367
bool m_do_qed_schwinger
Definition MultiParticleContainer.H:513
std::unique_ptr< amrex::Parser > m_By_particle_parser
Definition MultiParticleContainer.H:353
void PushX(amrex::Real dt)
This pushes the particle positions by one time step for all the species in the MultiParticleContainer...
Definition MultiParticleContainer.cpp:543
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_lengths
Definition MultiParticleContainer.H:364
std::string m_qed_schwinger_ele_product_name
Definition MultiParticleContainer.H:515
int m_nspecies_quantum_sync
Definition MultiParticleContainer.H:468
int m_qed_schwinger_ele_product
Definition MultiParticleContainer.H:519
int nSpecies() const
Definition MultiParticleContainer.H:283
int m_nspecies_breit_wheeler
Definition MultiParticleContainer.H:469
void MFItInfoCheckTiling(const WarpXParticleContainer &) const noexcept
Definition MultiParticleContainer.H:552
void InitQuantumSync()
Definition MultiParticleContainer.cpp:1279
std::shared_ptr< BreitWheelerEngine > m_shr_p_bw_engine
Definition MultiParticleContainer.H:457
void Redistribute()
Definition MultiParticleContainer.cpp:877
std::vector< PCTypes > species_types
Definition MultiParticleContainer.H:432
void PostRestart()
Definition MultiParticleContainer.cpp:450
int NSpeciesQuantumSync() const
Definition MultiParticleContainer.H:483
amrex::ParticleReal m_quantum_sync_photon_creation_energy_threshold
Definition MultiParticleContainer.H:477
void DepositCurrent(ablastr::fields::MultiLevelVectorField const &J, amrex::Real dt, amrex::Real relative_time)
Deposit current density.
Definition MultiParticleContainer.cpp:587
WarpXParticleContainer & GetParticleContainerFromName(const std::string &name) const
Definition MultiParticleContainer.cpp:404
std::vector< std::string > species_names
Definition MultiParticleContainer.H:419
amrex::Box ComputeSchwingerGlobalBox() const
Definition MultiParticleContainer.cpp:1692
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_strengths_B
Definition MultiParticleContainer.H:370
amrex::Real m_qed_schwinger_ymax
Definition MultiParticleContainer.H:535
void doQEDSchwinger()
Definition MultiParticleContainer.cpp:1584
std::string m_E_ext_particle_s
Definition MultiParticleContainer.H:350
amrex::Real m_qed_schwinger_xmin
Definition MultiParticleContainer.H:532
int nSpeciesGatherFromMainGrid() const
Definition MultiParticleContainer.H:313
WarpXParticleContainer * GetParticleContainerPtr(int index) const
Definition MultiParticleContainer.H:86
void defineAllParticleTiles()
Definition MultiParticleContainer.cpp:885
void BreitWheelerGenerateTable()
Definition MultiParticleContainer.cpp:1500
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:562
void ReadParameters()
Definition MultiParticleContainer.cpp:130
int doContinuousInjection() const
Definition MultiParticleContainer.cpp:1000
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_strengths_B
Definition MultiParticleContainer.H:366
int nSpeciesDepositOnMainGrid() const
Definition MultiParticleContainer.H:306
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:652
MultiParticleContainer & operator=(MultiParticleContainer &&)=default
std::unique_ptr< CollisionHandler > collisionhandler
Definition MultiParticleContainer.H:423
amrex::ParticleReal m_repeated_plasma_lens_period
Definition MultiParticleContainer.H:362
amrex::Real m_qed_schwinger_xmax
Definition MultiParticleContainer.H:533
void doCollisions(int step, amrex::Real cur_time, amrex::Real dt)
Definition MultiParticleContainer.cpp:1210
void MFItInfoCheckTiling(const WarpXParticleContainer &pc_src, First const &pc_dst, Args const &... others) const noexcept
Definition MultiParticleContainer.H:556
int getSpeciesID(const std::string &product_str) const
Definition MultiParticleContainer.cpp:1079
amrex::Vector< amrex::Long > NumberOfParticlesInGrid(int lev) const
Definition MultiParticleContainer.cpp:928
void InitBreitWheeler()
Definition MultiParticleContainer.cpp:1351
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:1144
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:919
std::unique_ptr< amrex::MultiFab > GetGlobalPlasmaFrequency(int lev)
Definition MultiParticleContainer.cpp:704
bool m_do_back_transformed_particles
Definition MultiParticleContainer.H:550
void GenerateGlobalDebyeLength()
Definition MultiParticleContainer.cpp:745
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, MomentumPushType momentum_push_type)
Definition MultiParticleContainer.cpp:551
void TransformMomentumToCurvilinear(bool forward)
Definition MultiParticleContainer.cpp:425
std::array< amrex::ParticleReal, 3 > meanParticleVelocity(int index)
Definition MultiParticleContainer.H:91
amrex::Vector< std::unique_ptr< WarpXParticleContainer > > allcontainers
Definition MultiParticleContainer.H:544
void SetParticleBoxArray(int lev, amrex::BoxArray &new_ba)
Definition MultiParticleContainer.cpp:959
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:1851
WarpXParticleContainer & GetParticleContainer(int index) const
Definition MultiParticleContainer.H:83
void ScrapeParticlesAtEB(ablastr::fields::MultiLevelScalarField const &distance_to_eb)
Definition MultiParticleContainer.cpp:1239
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:990
MultiParticleContainer(MultiParticleContainer &&)=default
void CheckIonizationProductSpecies()
Definition MultiParticleContainer.cpp:1228
std::unique_ptr< amrex::Parser > m_Ex_particle_parser
Definition MultiParticleContainer.H:356
std::vector< std::string > GetSpeciesAndLasersNames() const
Definition MultiParticleContainer.H:340
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_starts
Definition MultiParticleContainer.H:363
int m_qed_schwinger_threshold_poisson_gaussian
Definition MultiParticleContainer.H:528
std::vector< bool > m_gather_from_main_grid
instead of gathering fields from the finest patch level, gather from the coarsest
Definition MultiParticleContainer.H:430
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:1754
void DepositCharge(const ablastr::fields::MultiLevelScalarField &rho, amrex::Real relative_time)
Deposit charge density.
Definition MultiParticleContainer.cpp:615
void CheckQEDProductSpecies()
Definition MultiParticleContainer.cpp:1932
amrex::Real m_qed_schwinger_y_size
Definition MultiParticleContainer.H:523
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_strengths_E
Definition MultiParticleContainer.H:369
amrex::ParticleReal maxParticleVelocity()
Definition MultiParticleContainer.cpp:416
int nContainers() const
Definition MultiParticleContainer.H:285
ExternalParticleFields m_external_particle_fields_metadata
Definition MultiParticleContainer.H:360
void WriteHeader(std::ostream &os) const
Definition ParticleIO.cpp:242
std::string m_qed_schwinger_pos_product_name
Definition MultiParticleContainer.H:517
void RedistributeLocal(int max_cells_travelled)
Definition MultiParticleContainer.cpp:901
void Restart(const std::string &dir)
Definition ParticleIO.cpp:129
std::vector< std::string > GetSpeciesNames() const
Definition MultiParticleContainer.H:336
void ContinuousFluxInjection(amrex::Real t, amrex::Real dt) const
Definition MultiParticleContainer.cpp:1016
int m_qed_schwinger_pos_product
Definition MultiParticleContainer.H:521
amrex::MFItInfo getMFItInfo(const WarpXParticleContainer &pc_src, Args const &... pc_dsts) const noexcept
Definition MultiParticleContainer.H:436
MultiParticleContainer(amrex::AmrCore *amr_core)
Definition MultiParticleContainer.cpp:96
std::unique_ptr< amrex::Parser > m_Ez_particle_parser
Definition MultiParticleContainer.H:358
std::shared_ptr< QuantumSynchrotronEngine > m_shr_p_qs_engine
Definition MultiParticleContainer.H:458
~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:426
void SetDoBackTransformedParticles(bool do_back_transformed_particles)
Definition MultiParticleContainer.cpp:1103
PCTypes
Definition MultiParticleContainer.H:417
@ RigidInjected
Definition MultiParticleContainer.H:417
@ Physical
Definition MultiParticleContainer.H:417
@ Photon
Definition MultiParticleContainer.H:417
std::unique_ptr< amrex::Parser > m_Bx_particle_parser
Definition MultiParticleContainer.H:352
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:1768
Definition WarpXParticleContainer.H:195
amrex_real Real
amrex_particle_real ParticleReal
PODVector< T, ArenaAllocator< T > > DeviceVector
constexpr auto c2
square of the vacuum speed of light [m^2/s^2]
Definition constant.H:221
constexpr auto m_e
electron mass [kg]
Definition constant.H:172
amrex::Vector< ScalarField > MultiLevelScalarField
Definition MultiFabRegister.H:210
amrex::Vector< VectorField > MultiLevelVectorField
Definition MultiFabRegister.H:218
bool notInLaunchRegion() noexcept
BoxND< 3 > Box
IntVectND< 3 > IntVect
Definition ImplicitOptions.H:7
Definition MultiFabRegister.H:272
MFItInfo & SetDynamic(bool f) noexcept
MFItInfo & EnableTiling(const IntVect &ts=FabArrayBase::mfiter_tile_size) noexcept