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 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 int lev, amrex::Real dt);
126
131 void PushX (amrex::Real dt);
132
139 void PushP (int lev, amrex::Real dt,
140 const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
141 const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz,
142 MomentumPushType momentum_push_type);
143
150 std::unique_ptr<amrex::MultiFab> GetZeroChargeDensity(int lev);
151
161 void
163 amrex::Real relative_time);
164
176 void
178 amrex::Real dt, amrex::Real relative_time);
179
190 void
192
196
197 std::unique_ptr<amrex::MultiFab> GetChargeDensity(int lev, bool local = false);
198
200
201 void doFieldIonization (int lev,
202 const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
203 const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz);
204
205 void doCollisions (int step, amrex::Real cur_time, amrex::Real dt);
206
214 void doResampling (const amrex::Vector<amrex::Geometry>& geom, int timestep, bool verbose);
215
216#ifdef WARPX_QED
224 void doQEDSchwinger ();
225
230 [[nodiscard]] amrex::Box ComputeSchwingerGlobalBox () const;
231#endif
232
233 void Restart (const std::string& dir);
234
235 void PostRestart ();
236
237 void ReadHeader (std::istream& is);
238
239 void WriteHeader (std::ostream& os) const;
240
241 void SortParticlesByBin (
242 const amrex::IntVect& bin_size,
243 bool sort_particles_for_deposition,
244 const amrex::IntVect& sort_idx_type);
245
246 void Redistribute ();
247
249
251
252 void RedistributeLocal (int max_cells_travelled);
253
257
265 [[nodiscard]] amrex::Vector<amrex::Long> GetZeroParticlesInGrid(int lev) const;
266
267 [[nodiscard]] amrex::Vector<amrex::Long> NumberOfParticlesInGrid(int lev) const;
268
269 void Increment (amrex::MultiFab& mf, int lev);
270
271 void SetParticleBoxArray (int lev, amrex::BoxArray& new_ba);
273
274 [[nodiscard]] int nSpecies () const {return static_cast<int>(species_names.size());}
275 [[nodiscard]] int nLasers () const {return static_cast<int>(lasers_names.size());}
276 [[nodiscard]] int nContainers () const {return static_cast<int>(allcontainers.size());}
277
282 void SetDoBackTransformedParticles (bool do_back_transformed_particles);
288 void SetDoBackTransformedParticles (const std::string& species_name, bool do_back_transformed_particles);
289
292 [[nodiscard]] bool getDoBackTransformedParticles () const
293 {
295 }
296
297 [[nodiscard]] int nSpeciesDepositOnMainGrid () const
298 {
299 bool const onMainGrid = true;
300 auto const & v = m_deposit_on_main_grid;
301 return static_cast<int>(std::count( v.begin(), v.end(), onMainGrid ));
302 }
303
304 [[nodiscard]] int nSpeciesGatherFromMainGrid() const
305 {
306 bool const fromMainGrid = true;
307 auto const & v = m_gather_from_main_grid;
308 return static_cast<int>(std::count( v.begin(), v.end(), fromMainGrid ));
309 }
310
311 // Inject particles during the simulation (for particles entering the
312 // simulation domain after some iterations, due to flowing plasma and/or
313 // moving window).
314 void ContinuousInjection(const amrex::RealBox& injection_box) const;
315
320 void UpdateAntennaPosition(amrex::Real dt) const;
321
322 [[nodiscard]] int doContinuousInjection() const;
323
324 // Inject particles from a surface during the simulation
326
327 [[nodiscard]] std::vector<std::string> GetSpeciesNames() const { return species_names; }
328
329 [[nodiscard]] std::vector<std::string> GetLasersNames() const { return lasers_names; }
330
331 [[nodiscard]] std::vector<std::string> GetSpeciesAndLasersNames() const
332 {
333 std::vector<std::string> tmp = species_names;
334 tmp.insert(tmp.end(), lasers_names.begin(), lasers_names.end());
335 return tmp;
336 }
337
339
340 std::string m_B_ext_particle_s = "none";
341 std::string m_E_ext_particle_s = "none";
342 // Parser for B_external on the particle
343 std::unique_ptr<amrex::Parser> m_Bx_particle_parser;
344 std::unique_ptr<amrex::Parser> m_By_particle_parser;
345 std::unique_ptr<amrex::Parser> m_Bz_particle_parser;
346 // Parser for E_external on the particle
347 std::unique_ptr<amrex::Parser> m_Ex_particle_parser;
348 std::unique_ptr<amrex::Parser> m_Ey_particle_parser;
349 std::unique_ptr<amrex::Parser> m_Ez_particle_parser;
350
352
362
363#ifdef WARPX_QED
367 void doQedEvents (int lev,
368 const amrex::MultiFab& Ex,
369 const amrex::MultiFab& Ey,
370 const amrex::MultiFab& Ez,
371 const amrex::MultiFab& Bx,
372 const amrex::MultiFab& By,
373 const amrex::MultiFab& Bz);
374#endif
375
376 [[nodiscard]] int getSpeciesID (const std::string& product_str) const;
377
380
381protected:
382
383#ifdef WARPX_QED
387 void doQedBreitWheeler (int lev,
388 const amrex::MultiFab& Ex,
389 const amrex::MultiFab& Ey,
390 const amrex::MultiFab& Ez,
391 const amrex::MultiFab& Bx,
392 const amrex::MultiFab& By,
393 const amrex::MultiFab& Bz);
394
398 void doQedQuantumSync (int lev,
399 const amrex::MultiFab& Ex,
400 const amrex::MultiFab& Ey,
401 const amrex::MultiFab& Ez,
402 const amrex::MultiFab& Bx,
403 const amrex::MultiFab& By,
404 const amrex::MultiFab& Bz);
405#endif
406
407 // Particle container types
409
410 std::vector<std::string> species_names;
411
412 std::vector<std::string> lasers_names;
413
414 std::unique_ptr<CollisionHandler> collisionhandler;
415
417 std::vector<bool> m_deposit_on_main_grid;
419
421 std::vector<bool> m_gather_from_main_grid;
422
423 std::vector<PCTypes> species_types;
424
425 template<typename ...Args>
426 [[nodiscard]]
428 Args const&... pc_dsts) const noexcept
429 {
430 amrex::MFItInfo info;
431
432 MFItInfoCheckTiling(pc_src, pc_dsts...);
433
436 }
437
438#ifdef AMREX_USE_OMP
439 info.SetDynamic(true);
440#endif
441
442 return info;
443 }
444
445
446#ifdef WARPX_QED
447 // The QED engines
448 std::shared_ptr<BreitWheelerEngine> m_shr_p_bw_engine;
449 std::shared_ptr<QuantumSynchrotronEngine> m_shr_p_qs_engine;
450 //_______________________________
451
456 void InitQED ();
457
458 //Variables to store how many species need a QED process
461 //________
462
464 static_cast<amrex::ParticleReal>(
465 2.0 * PhysConst::m_e * PhysConst::c2 );
466
467
470
474 [[nodiscard]] int NSpeciesQuantumSync() const { return m_nspecies_quantum_sync;}
475
479 [[nodiscard]] int NSpeciesBreitWheeler() const { return m_nspecies_breit_wheeler;}
480
484 void InitQuantumSync ();
485
489 void InitBreitWheeler ();
490
496
502
504 bool m_do_qed_schwinger = false;
523 amrex::Real m_qed_schwinger_xmin = std::numeric_limits<amrex::Real>::lowest();
524 amrex::Real m_qed_schwinger_xmax = std::numeric_limits<amrex::Real>::max();
525 amrex::Real m_qed_schwinger_ymin = std::numeric_limits<amrex::Real>::lowest();
526 amrex::Real m_qed_schwinger_ymax = std::numeric_limits<amrex::Real>::max();
527 amrex::Real m_qed_schwinger_zmin = std::numeric_limits<amrex::Real>::lowest();
528 amrex::Real m_qed_schwinger_zmax = std::numeric_limits<amrex::Real>::max();
529
530#endif
531
532private:
533
534 // physical particles (+ laser)
536
537 void ReadParameters ();
538
539 void mapSpeciesProduct ();
540
542
543 void MFItInfoCheckTiling(const WarpXParticleContainer& /*pc_src*/) const noexcept
544 {}
545
546 template<typename First, typename ...Args>
548 First const& pc_dst, Args const&... others) const noexcept
549 {
551 WARPX_ALWAYS_ASSERT_WITH_MESSAGE(pc_dst.do_tiling,
552 "For particle creation processes, either all or none of the "
553 "particle species must use tiling.");
554 }
555
556 MFItInfoCheckTiling(pc_src, others...);
557 }
558
565
566#ifdef WARPX_QED
573#endif
574
575
576};
577#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:772
std::vector< std::string > GetLasersNames() const
Definition MultiParticleContainer.H:329
amrex::Vector< std::unique_ptr< WarpXParticleContainer > >::iterator begin()
Definition MultiParticleContainer.H:378
void deleteInvalidParticles()
Definition MultiParticleContainer.cpp:804
std::unique_ptr< amrex::Parser > m_Bz_particle_parser
Definition MultiParticleContainer.H:345
void InitQED()
Definition MultiParticleContainer.cpp:1159
std::unique_ptr< amrex::Parser > m_Ey_particle_parser
Definition MultiParticleContainer.H:348
std::unique_ptr< amrex::MultiFab > GetChargeDensity(int lev, bool local=false)
Definition MultiParticleContainer.cpp:675
int nLasers() const
Definition MultiParticleContainer.H:275
bool getDoBackTransformedParticles() const
Definition MultiParticleContainer.H:292
amrex::Real m_qed_schwinger_zmax
Definition MultiParticleContainer.H:528
amrex::Real m_qed_schwinger_zmin
Definition MultiParticleContainer.H:527
void InitMultiPhysicsModules()
Definition MultiParticleContainer.cpp:453
std::string m_B_ext_particle_s
Definition MultiParticleContainer.H:340
void QuantumSyncGenerateTable()
Definition MultiParticleContainer.cpp:1322
MultiParticleContainer & operator=(MultiParticleContainer const &)=delete
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_strengths_E
Definition MultiParticleContainer.H:356
void AllocData()
Definition MultiParticleContainer.cpp:425
void SetParticleDistributionMap(int lev, amrex::DistributionMapping &new_dm)
Definition MultiParticleContainer.cpp:878
void InitData()
Definition MultiParticleContainer.cpp:433
amrex::Real m_qed_schwinger_ymin
Definition MultiParticleContainer.H:525
std::vector< bool > m_laser_deposit_on_main_grid
Definition MultiParticleContainer.H:418
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:471
void Increment(amrex::MultiFab &mf, int lev)
Definition MultiParticleContainer.cpp:862
void mapSpeciesProduct()
Definition MultiParticleContainer.cpp:939
void ContinuousInjection(const amrex::RealBox &injection_box) const
Definition MultiParticleContainer.cpp:891
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_lengths
Definition MultiParticleContainer.H:359
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:1127
void DepositMassMatrices(ablastr::fields::MultiFabRegister &fields, int lev, amrex::Real dt)
Deposit mass matrices.
Definition MultiParticleContainer.cpp:519
amrex::Vector< std::unique_ptr< WarpXParticleContainer > >::iterator end()
Definition MultiParticleContainer.H:379
int NSpeciesBreitWheeler() const
Definition MultiParticleContainer.H:479
void ReadHeader(std::istream &is)
Definition ParticleIO.cpp:231
static constexpr auto m_default_quantum_sync_photon_creation_energy_threshold
Definition MultiParticleContainer.H:463
void ApplyBoundaryConditions()
Definition MultiParticleContainer.cpp:822
std::vector< std::string > lasers_names
Definition MultiParticleContainer.H:412
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_starts
Definition MultiParticleContainer.H:358
bool m_do_qed_schwinger
Definition MultiParticleContainer.H:504
std::unique_ptr< amrex::Parser > m_By_particle_parser
Definition MultiParticleContainer.H:344
void PushX(amrex::Real dt)
This pushes the particle positions by one time step for all the species in the MultiParticleContainer...
Definition MultiParticleContainer.cpp:536
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_lengths
Definition MultiParticleContainer.H:355
std::string m_qed_schwinger_ele_product_name
Definition MultiParticleContainer.H:506
int m_nspecies_quantum_sync
Definition MultiParticleContainer.H:459
int m_qed_schwinger_ele_product
Definition MultiParticleContainer.H:510
int nSpecies() const
Definition MultiParticleContainer.H:274
int m_nspecies_breit_wheeler
Definition MultiParticleContainer.H:460
void MFItInfoCheckTiling(const WarpXParticleContainer &) const noexcept
Definition MultiParticleContainer.H:543
void InitQuantumSync()
Definition MultiParticleContainer.cpp:1190
std::shared_ptr< BreitWheelerEngine > m_shr_p_bw_engine
Definition MultiParticleContainer.H:448
void Redistribute()
Definition MultiParticleContainer.cpp:788
std::vector< PCTypes > species_types
Definition MultiParticleContainer.H:423
void PostRestart()
Definition MultiParticleContainer.cpp:443
int NSpeciesQuantumSync() const
Definition MultiParticleContainer.H:474
amrex::ParticleReal m_quantum_sync_photon_creation_energy_threshold
Definition MultiParticleContainer.H:468
void DepositCurrent(ablastr::fields::MultiLevelVectorField const &J, amrex::Real dt, amrex::Real relative_time)
Deposit current density.
Definition MultiParticleContainer.cpp:580
WarpXParticleContainer & GetParticleContainerFromName(const std::string &name) const
Definition MultiParticleContainer.cpp:404
std::vector< std::string > species_names
Definition MultiParticleContainer.H:410
amrex::Box ComputeSchwingerGlobalBox() const
Definition MultiParticleContainer.cpp:1603
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_strengths_B
Definition MultiParticleContainer.H:361
amrex::Real m_qed_schwinger_ymax
Definition MultiParticleContainer.H:526
void doQEDSchwinger()
Definition MultiParticleContainer.cpp:1495
std::string m_E_ext_particle_s
Definition MultiParticleContainer.H:341
amrex::Real m_qed_schwinger_xmin
Definition MultiParticleContainer.H:523
int nSpeciesGatherFromMainGrid() const
Definition MultiParticleContainer.H:304
WarpXParticleContainer * GetParticleContainerPtr(int index) const
Definition MultiParticleContainer.H:86
void defineAllParticleTiles()
Definition MultiParticleContainer.cpp:796
void BreitWheelerGenerateTable()
Definition MultiParticleContainer.cpp:1411
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:555
void ReadParameters()
Definition MultiParticleContainer.cpp:130
int doContinuousInjection() const
Definition MultiParticleContainer.cpp:911
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_strengths_B
Definition MultiParticleContainer.H:357
int nSpeciesDepositOnMainGrid() const
Definition MultiParticleContainer.H:297
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:645
MultiParticleContainer & operator=(MultiParticleContainer &&)=default
std::unique_ptr< CollisionHandler > collisionhandler
Definition MultiParticleContainer.H:414
amrex::ParticleReal m_repeated_plasma_lens_period
Definition MultiParticleContainer.H:353
amrex::Real m_qed_schwinger_xmax
Definition MultiParticleContainer.H:524
void doCollisions(int step, amrex::Real cur_time, amrex::Real dt)
Definition MultiParticleContainer.cpp:1121
void MFItInfoCheckTiling(const WarpXParticleContainer &pc_src, First const &pc_dst, Args const &... others) const noexcept
Definition MultiParticleContainer.H:547
int getSpeciesID(const std::string &product_str) const
Definition MultiParticleContainer.cpp:990
amrex::Vector< amrex::Long > NumberOfParticlesInGrid(int lev) const
Definition MultiParticleContainer.cpp:839
void InitBreitWheeler()
Definition MultiParticleContainer.cpp:1262
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:1055
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:830
bool m_do_back_transformed_particles
Definition MultiParticleContainer.H:541
void GenerateGlobalDebyeLength()
Definition MultiParticleContainer.cpp:697
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:544
std::array< amrex::ParticleReal, 3 > meanParticleVelocity(int index)
Definition MultiParticleContainer.H:91
amrex::Vector< std::unique_ptr< WarpXParticleContainer > > allcontainers
Definition MultiParticleContainer.H:535
void SetParticleBoxArray(int lev, amrex::BoxArray &new_ba)
Definition MultiParticleContainer.cpp:870
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:1762
WarpXParticleContainer & GetParticleContainer(int index) const
Definition MultiParticleContainer.H:83
void ScrapeParticlesAtEB(ablastr::fields::MultiLevelScalarField const &distance_to_eb)
Definition MultiParticleContainer.cpp:1150
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:901
MultiParticleContainer(MultiParticleContainer &&)=default
void CheckIonizationProductSpecies()
Definition MultiParticleContainer.cpp:1139
std::unique_ptr< amrex::Parser > m_Ex_particle_parser
Definition MultiParticleContainer.H:347
std::vector< std::string > GetSpeciesAndLasersNames() const
Definition MultiParticleContainer.H:331
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_starts
Definition MultiParticleContainer.H:354
int m_qed_schwinger_threshold_poisson_gaussian
Definition MultiParticleContainer.H:519
std::vector< bool > m_gather_from_main_grid
instead of gathering fields from the finest patch level, gather from the coarsest
Definition MultiParticleContainer.H:421
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:1665
void DepositCharge(const ablastr::fields::MultiLevelScalarField &rho, amrex::Real relative_time)
Deposit charge density.
Definition MultiParticleContainer.cpp:608
void CheckQEDProductSpecies()
Definition MultiParticleContainer.cpp:1843
amrex::Real m_qed_schwinger_y_size
Definition MultiParticleContainer.H:514
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_strengths_E
Definition MultiParticleContainer.H:360
amrex::ParticleReal maxParticleVelocity()
Definition MultiParticleContainer.cpp:416
int nContainers() const
Definition MultiParticleContainer.H:276
ExternalParticleFields m_external_particle_fields_metadata
Definition MultiParticleContainer.H:351
void WriteHeader(std::ostream &os) const
Definition ParticleIO.cpp:242
std::string m_qed_schwinger_pos_product_name
Definition MultiParticleContainer.H:508
void RedistributeLocal(int max_cells_travelled)
Definition MultiParticleContainer.cpp:812
void Restart(const std::string &dir)
Definition ParticleIO.cpp:129
std::vector< std::string > GetSpeciesNames() const
Definition MultiParticleContainer.H:327
void ContinuousFluxInjection(amrex::Real t, amrex::Real dt) const
Definition MultiParticleContainer.cpp:927
int m_qed_schwinger_pos_product
Definition MultiParticleContainer.H:512
amrex::MFItInfo getMFItInfo(const WarpXParticleContainer &pc_src, Args const &... pc_dsts) const noexcept
Definition MultiParticleContainer.H:427
MultiParticleContainer(amrex::AmrCore *amr_core)
Definition MultiParticleContainer.cpp:96
std::unique_ptr< amrex::Parser > m_Ez_particle_parser
Definition MultiParticleContainer.H:349
std::shared_ptr< QuantumSynchrotronEngine > m_shr_p_qs_engine
Definition MultiParticleContainer.H:449
~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:417
void SetDoBackTransformedParticles(bool do_back_transformed_particles)
Definition MultiParticleContainer.cpp:1014
PCTypes
Definition MultiParticleContainer.H:408
@ RigidInjected
Definition MultiParticleContainer.H:408
@ Physical
Definition MultiParticleContainer.H:408
@ Photon
Definition MultiParticleContainer.H:408
std::unique_ptr< amrex::Parser > m_Bx_particle_parser
Definition MultiParticleContainer.H:343
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:1679
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:214
constexpr auto m_e
electron mass [kg]
Definition constant.H:165
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