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"
27
29
30#include <AMReX_BLassert.H>
31#include <AMReX_Box.H>
32#include <AMReX_Config.H>
33#include <AMReX_Geometry.H>
34#include <AMReX_GpuControl.H>
35#include <AMReX_INT.H>
36#include <AMReX_MFIter.H>
37#include <AMReX_REAL.H>
38#include <AMReX_RealBox.H>
39#include <AMReX_Vector.H>
40
41#include <AMReX_BaseFwd.H>
42#include <AMReX_AmrCoreFwd.H>
43
44#include <algorithm>
45#include <array>
46#include <iosfwd>
47#include <iterator>
48#include <limits>
49#include <memory>
50#include <string>
51#include <vector>
52
68{
69
70public:
71
73
75
80
81 [[nodiscard]] WarpXParticleContainer&
82 GetParticleContainer (int index) const {return *allcontainers[index];}
83
84 [[nodiscard]] WarpXParticleContainer*
85 GetParticleContainerPtr (int index) const {return allcontainers[index].get();}
86
87 [[nodiscard]] WarpXParticleContainer&
88 GetParticleContainerFromName (const std::string& name) const;
89
90 std::array<amrex::ParticleReal, 3> meanParticleVelocity(int index) {
91 return allcontainers[index]->meanParticleVelocity();
92 }
93
94 amrex::ParticleReal maxParticleVelocity();
95
96 void AllocData ();
97
98 void InitData ();
99
101
107 void Evolve (
109 int lev,
110 std::string const& current_fp_string,
111 amrex::Real t,
112 amrex::Real dt,
113 SubcyclingHalf subcycling_half=SubcyclingHalf::None,
114 bool skip_deposition=false,
115 ImplicitOptions const * implicit_options = nullptr
116 );
117
122 void PushX (amrex::Real dt);
123
130 void PushP (int lev, amrex::Real dt,
131 const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
132 const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz);
133
140 std::unique_ptr<amrex::MultiFab> GetZeroChargeDensity(int lev);
141
151 void
153 amrex::Real relative_time);
154
166 void
168 amrex::Real dt, amrex::Real relative_time);
169
180 void
181 DepositTemperatures (ablastr::fields::MultiFabRegister & fields, amrex::Real relative_time);
182
186
187 std::unique_ptr<amrex::MultiFab> GetChargeDensity(int lev, bool local = false);
188
190
191 void doFieldIonization (int lev,
192 const amrex::MultiFab& Ex, const amrex::MultiFab& Ey, const amrex::MultiFab& Ez,
193 const amrex::MultiFab& Bx, const amrex::MultiFab& By, const amrex::MultiFab& Bz);
194
195 void doCollisions (int step, amrex::Real cur_time, amrex::Real dt);
196
204 void doResampling (const amrex::Vector<amrex::Geometry>& geom, int timestep, bool verbose);
205
206#ifdef WARPX_QED
214 void doQEDSchwinger ();
215
220 [[nodiscard]] amrex::Box ComputeSchwingerGlobalBox () const;
221#endif
222
223 void Restart (const std::string& dir);
224
225 void PostRestart ();
226
227 void ReadHeader (std::istream& is);
228
229 void WriteHeader (std::ostream& os) const;
230
231 void SortParticlesByBin (
232 const amrex::IntVect& bin_size,
233 bool sort_particles_for_deposition,
234 const amrex::IntVect& sort_idx_type);
235
236 void Redistribute ();
237
239
241
242 void RedistributeLocal (int num_ghost);
243
247
255 [[nodiscard]] amrex::Vector<amrex::Long> GetZeroParticlesInGrid(int lev) const;
256
257 [[nodiscard]] amrex::Vector<amrex::Long> NumberOfParticlesInGrid(int lev) const;
258
259 void Increment (amrex::MultiFab& mf, int lev);
260
261 void SetParticleBoxArray (int lev, amrex::BoxArray& new_ba);
263
264 [[nodiscard]] int nSpecies () const {return static_cast<int>(species_names.size());}
265 [[nodiscard]] int nLasers () const {return static_cast<int>(lasers_names.size());}
266 [[nodiscard]] int nContainers () const {return static_cast<int>(allcontainers.size());}
267
272 void SetDoBackTransformedParticles (bool do_back_transformed_particles);
278 void SetDoBackTransformedParticles (const std::string& species_name, bool do_back_transformed_particles);
279
282 [[nodiscard]] bool getDoBackTransformedParticles () const
283 {
285 }
286
287 [[nodiscard]] int nSpeciesDepositOnMainGrid () const
288 {
289 bool const onMainGrid = true;
290 auto const & v = m_deposit_on_main_grid;
291 return static_cast<int>(std::count( v.begin(), v.end(), onMainGrid ));
292 }
293
294 [[nodiscard]] int nSpeciesGatherFromMainGrid() const
295 {
296 bool const fromMainGrid = true;
297 auto const & v = m_gather_from_main_grid;
298 return static_cast<int>(std::count( v.begin(), v.end(), fromMainGrid ));
299 }
300
301 // Inject particles during the simulation (for particles entering the
302 // simulation domain after some iterations, due to flowing plasma and/or
303 // moving window).
304 void ContinuousInjection(const amrex::RealBox& injection_box) const;
305
310 void UpdateAntennaPosition(amrex::Real dt) const;
311
312 [[nodiscard]] int doContinuousInjection() const;
313
314 // Inject particles from a surface during the simulation
315 void ContinuousFluxInjection(amrex::Real t, amrex::Real dt) const;
316
317 [[nodiscard]] std::vector<std::string> GetSpeciesNames() const { return species_names; }
318
319 [[nodiscard]] std::vector<std::string> GetLasersNames() const { return lasers_names; }
320
321 [[nodiscard]] std::vector<std::string> GetSpeciesAndLasersNames() const
322 {
323 std::vector<std::string> tmp = species_names;
324 tmp.insert(tmp.end(), lasers_names.begin(), lasers_names.end());
325 return tmp;
326 }
327
329
330 std::string m_B_ext_particle_s = "none";
331 std::string m_E_ext_particle_s = "none";
332 // Parser for B_external on the particle
333 std::unique_ptr<amrex::Parser> m_Bx_particle_parser;
334 std::unique_ptr<amrex::Parser> m_By_particle_parser;
335 std::unique_ptr<amrex::Parser> m_Bz_particle_parser;
336 // Parser for E_external on the particle
337 std::unique_ptr<amrex::Parser> m_Ex_particle_parser;
338 std::unique_ptr<amrex::Parser> m_Ey_particle_parser;
339 std::unique_ptr<amrex::Parser> m_Ez_particle_parser;
340 // Parser for time dependency of B_external from file on the particle
341 std::unique_ptr<amrex::Parser> m_B_particle_from_file_parser;
342 // Parser for time dependency of B_external from file on the particle
343 std::unique_ptr<amrex::Parser> m_E_particle_from_file_parser;
344
347
348
349 amrex::ParticleReal m_repeated_plasma_lens_period;
358
359#ifdef WARPX_QED
363 void doQedEvents (int lev,
364 const amrex::MultiFab& Ex,
365 const amrex::MultiFab& Ey,
366 const amrex::MultiFab& Ez,
367 const amrex::MultiFab& Bx,
368 const amrex::MultiFab& By,
369 const amrex::MultiFab& Bz);
370#endif
371
372 [[nodiscard]] int getSpeciesID (const std::string& product_str) const;
373
376
377protected:
378
379#ifdef WARPX_QED
383 void doQedBreitWheeler (int lev,
384 const amrex::MultiFab& Ex,
385 const amrex::MultiFab& Ey,
386 const amrex::MultiFab& Ez,
387 const amrex::MultiFab& Bx,
388 const amrex::MultiFab& By,
389 const amrex::MultiFab& Bz);
390
394 void doQedQuantumSync (int lev,
395 const amrex::MultiFab& Ex,
396 const amrex::MultiFab& Ey,
397 const amrex::MultiFab& Ez,
398 const amrex::MultiFab& Bx,
399 const amrex::MultiFab& By,
400 const amrex::MultiFab& Bz);
401#endif
402
403 // Particle container types
405
406 std::vector<std::string> species_names;
407
408 std::vector<std::string> lasers_names;
409
410 std::unique_ptr<CollisionHandler> collisionhandler;
411
413 std::vector<bool> m_deposit_on_main_grid;
415
417 std::vector<bool> m_gather_from_main_grid;
418
419 std::vector<PCTypes> species_types;
420
421 template<typename ...Args>
422 [[nodiscard]]
424 Args const&... pc_dsts) const noexcept
425 {
426 amrex::MFItInfo info;
427
428 MFItInfoCheckTiling(pc_src, pc_dsts...);
429
432 }
433
434#ifdef AMREX_USE_OMP
435 info.SetDynamic(true);
436#endif
437
438 return info;
439 }
440
441
442#ifdef WARPX_QED
443 // The QED engines
444 std::shared_ptr<BreitWheelerEngine> m_shr_p_bw_engine;
445 std::shared_ptr<QuantumSynchrotronEngine> m_shr_p_qs_engine;
446 //_______________________________
447
452 void InitQED ();
453
454 //Variables to store how many species need a QED process
457 //________
458
460 static_cast<amrex::ParticleReal>(
462
463
466
470 [[nodiscard]] int NSpeciesQuantumSync() const { return m_nspecies_quantum_sync;}
471
475 [[nodiscard]] int NSpeciesBreitWheeler() const { return m_nspecies_breit_wheeler;}
476
480 void InitQuantumSync ();
481
485 void InitBreitWheeler ();
486
492
498
500 bool m_do_qed_schwinger = false;
519 amrex::Real m_qed_schwinger_xmin = std::numeric_limits<amrex::Real>::lowest();
520 amrex::Real m_qed_schwinger_xmax = std::numeric_limits<amrex::Real>::max();
521 amrex::Real m_qed_schwinger_ymin = std::numeric_limits<amrex::Real>::lowest();
522 amrex::Real m_qed_schwinger_ymax = std::numeric_limits<amrex::Real>::max();
523 amrex::Real m_qed_schwinger_zmin = std::numeric_limits<amrex::Real>::lowest();
524 amrex::Real m_qed_schwinger_zmax = std::numeric_limits<amrex::Real>::max();
525
526#endif
527
528private:
529
530 // physical particles (+ laser)
532
533 void ReadParameters ();
534
535 void mapSpeciesProduct ();
536
538
539 void MFItInfoCheckTiling(const WarpXParticleContainer& /*pc_src*/) const noexcept
540 {}
541
542 template<typename First, typename ...Args>
544 First const& pc_dst, Args const&... others) const noexcept
545 {
547 WARPX_ALWAYS_ASSERT_WITH_MESSAGE(pc_dst.do_tiling,
548 "For particle creation processes, either all or none of the "
549 "particle species must use tiling.");
550 }
551
552 MFItInfoCheckTiling(pc_src, others...);
553 }
554
561
562#ifdef WARPX_QED
569#endif
570
571
572};
573#endif /*WARPX_ParticleContainer_H_*/
@ v
Definition RigidInjectedParticleContainer.H:27
#define WARPX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition TextMsg.H:13
SubcyclingHalf
Subcycling half selector.
Definition WarpXAlgorithmSelection.H:166
@ None
Definition WarpXAlgorithmSelection.H:166
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:319
void RedistributeLocal(int num_ghost)
Definition MultiParticleContainer.cpp:812
amrex::Vector< std::unique_ptr< WarpXParticleContainer > >::iterator begin()
Definition MultiParticleContainer.H:374
void deleteInvalidParticles()
Definition MultiParticleContainer.cpp:804
std::unique_ptr< amrex::Parser > m_Bz_particle_parser
Definition MultiParticleContainer.H:335
void InitQED()
Definition MultiParticleContainer.cpp:1151
std::unique_ptr< amrex::Parser > m_Ey_particle_parser
Definition MultiParticleContainer.H:338
std::unique_ptr< amrex::MultiFab > GetChargeDensity(int lev, bool local=false)
Definition MultiParticleContainer.cpp:672
int nLasers() const
Definition MultiParticleContainer.H:265
bool getDoBackTransformedParticles() const
Definition MultiParticleContainer.H:282
amrex::Real m_qed_schwinger_zmax
Definition MultiParticleContainer.H:524
amrex::Real m_qed_schwinger_zmin
Definition MultiParticleContainer.H:523
void InitMultiPhysicsModules()
Definition MultiParticleContainer.cpp:473
std::string m_B_ext_particle_s
Definition MultiParticleContainer.H:330
void QuantumSyncGenerateTable()
Definition MultiParticleContainer.cpp:1314
MultiParticleContainer & operator=(MultiParticleContainer const &)=delete
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_strengths_E
Definition MultiParticleContainer.H:352
void AllocData()
Definition MultiParticleContainer.cpp:445
void SetParticleDistributionMap(int lev, amrex::DistributionMapping &new_dm)
Definition MultiParticleContainer.cpp:876
void InitData()
Definition MultiParticleContainer.cpp:453
amrex::Real m_qed_schwinger_ymin
Definition MultiParticleContainer.H:521
std::vector< bool > m_laser_deposit_on_main_grid
Definition MultiParticleContainer.H:414
void Increment(amrex::MultiFab &mf, int lev)
Definition MultiParticleContainer.cpp:860
void mapSpeciesProduct()
Definition MultiParticleContainer.cpp:937
void ContinuousInjection(const amrex::RealBox &injection_box) const
Definition MultiParticleContainer.cpp:889
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_lengths
Definition MultiParticleContainer.H:355
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:1119
amrex::Vector< std::unique_ptr< WarpXParticleContainer > >::iterator end()
Definition MultiParticleContainer.H:375
int NSpeciesBreitWheeler() const
Definition MultiParticleContainer.H:475
void ReadHeader(std::istream &is)
Definition ParticleIO.cpp:231
static constexpr auto m_default_quantum_sync_photon_creation_energy_threshold
Definition MultiParticleContainer.H:459
void ApplyBoundaryConditions()
Definition MultiParticleContainer.cpp:820
std::vector< std::string > lasers_names
Definition MultiParticleContainer.H:408
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_starts
Definition MultiParticleContainer.H:354
bool m_do_qed_schwinger
Definition MultiParticleContainer.H:500
std::unique_ptr< amrex::Parser > m_By_particle_parser
Definition MultiParticleContainer.H:334
void PushX(amrex::Real dt)
This pushes the particle positions by one time step for all the species in the MultiParticleContainer...
Definition MultiParticleContainer.cpp:534
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_lengths
Definition MultiParticleContainer.H:351
std::string m_qed_schwinger_ele_product_name
Definition MultiParticleContainer.H:502
int m_nspecies_quantum_sync
Definition MultiParticleContainer.H:455
int m_qed_schwinger_ele_product
Definition MultiParticleContainer.H:506
int nSpecies() const
Definition MultiParticleContainer.H:264
int m_nspecies_breit_wheeler
Definition MultiParticleContainer.H:456
void MFItInfoCheckTiling(const WarpXParticleContainer &) const noexcept
Definition MultiParticleContainer.H:539
void InitQuantumSync()
Definition MultiParticleContainer.cpp:1182
std::shared_ptr< BreitWheelerEngine > m_shr_p_bw_engine
Definition MultiParticleContainer.H:444
void Redistribute()
Definition MultiParticleContainer.cpp:788
std::vector< PCTypes > species_types
Definition MultiParticleContainer.H:419
void PostRestart()
Definition MultiParticleContainer.cpp:463
int NSpeciesQuantumSync() const
Definition MultiParticleContainer.H:470
amrex::ParticleReal m_quantum_sync_photon_creation_energy_threshold
Definition MultiParticleContainer.H:464
void DepositCurrent(ablastr::fields::MultiLevelVectorField const &J, amrex::Real dt, amrex::Real relative_time)
Deposit current density.
Definition MultiParticleContainer.cpp:577
WarpXParticleContainer & GetParticleContainerFromName(const std::string &name) const
Definition MultiParticleContainer.cpp:425
std::vector< std::string > species_names
Definition MultiParticleContainer.H:406
amrex::Box ComputeSchwingerGlobalBox() const
Definition MultiParticleContainer.cpp:1597
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_strengths_B
Definition MultiParticleContainer.H:357
amrex::Real m_qed_schwinger_ymax
Definition MultiParticleContainer.H:522
void doQEDSchwinger()
Definition MultiParticleContainer.cpp:1489
std::string m_E_ext_particle_s
Definition MultiParticleContainer.H:331
amrex::Real m_qed_schwinger_xmin
Definition MultiParticleContainer.H:519
int nSpeciesGatherFromMainGrid() const
Definition MultiParticleContainer.H:294
amrex::ParserExecutor< 1 > m_Bfield_time_partparser
Definition MultiParticleContainer.H:346
WarpXParticleContainer * GetParticleContainerPtr(int index) const
Definition MultiParticleContainer.H:85
void defineAllParticleTiles()
Definition MultiParticleContainer.cpp:796
void BreitWheelerGenerateTable()
Definition MultiParticleContainer.cpp:1404
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:552
void ReadParameters()
Definition MultiParticleContainer.cpp:128
int doContinuousInjection() const
Definition MultiParticleContainer.cpp:909
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_strengths_B
Definition MultiParticleContainer.H:353
int nSpeciesDepositOnMainGrid() const
Definition MultiParticleContainer.H:287
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:642
MultiParticleContainer & operator=(MultiParticleContainer &&)=default
std::unique_ptr< CollisionHandler > collisionhandler
Definition MultiParticleContainer.H:410
amrex::ParticleReal m_repeated_plasma_lens_period
Definition MultiParticleContainer.H:349
amrex::Real m_qed_schwinger_xmax
Definition MultiParticleContainer.H:520
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:542
void doCollisions(int step, amrex::Real cur_time, amrex::Real dt)
Definition MultiParticleContainer.cpp:1113
void MFItInfoCheckTiling(const WarpXParticleContainer &pc_src, First const &pc_dst, Args const &... others) const noexcept
Definition MultiParticleContainer.H:543
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, ImplicitOptions const *implicit_options=nullptr)
This evolves all the particles by one PIC time step, including current deposition,...
Definition MultiParticleContainer.cpp:491
int getSpeciesID(const std::string &product_str) const
Definition MultiParticleContainer.cpp:982
amrex::Vector< amrex::Long > NumberOfParticlesInGrid(int lev) const
Definition MultiParticleContainer.cpp:837
void InitBreitWheeler()
Definition MultiParticleContainer.cpp:1254
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:1047
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:828
bool m_do_back_transformed_particles
Definition MultiParticleContainer.H:537
void GenerateGlobalDebyeLength()
Definition MultiParticleContainer.cpp:694
std::array< amrex::ParticleReal, 3 > meanParticleVelocity(int index)
Definition MultiParticleContainer.H:90
amrex::Vector< std::unique_ptr< WarpXParticleContainer > > allcontainers
Definition MultiParticleContainer.H:531
void SetParticleBoxArray(int lev, amrex::BoxArray &new_ba)
Definition MultiParticleContainer.cpp:868
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:1756
WarpXParticleContainer & GetParticleContainer(int index) const
Definition MultiParticleContainer.H:82
void ScrapeParticlesAtEB(ablastr::fields::MultiLevelScalarField const &distance_to_eb)
Definition MultiParticleContainer.cpp:1142
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:899
MultiParticleContainer(MultiParticleContainer &&)=default
void CheckIonizationProductSpecies()
Definition MultiParticleContainer.cpp:1131
std::unique_ptr< amrex::Parser > m_Ex_particle_parser
Definition MultiParticleContainer.H:337
std::vector< std::string > GetSpeciesAndLasersNames() const
Definition MultiParticleContainer.H:321
amrex::Vector< amrex::ParticleReal > h_repeated_plasma_lens_starts
Definition MultiParticleContainer.H:350
int m_qed_schwinger_threshold_poisson_gaussian
Definition MultiParticleContainer.H:515
std::vector< bool > m_gather_from_main_grid
instead of gathering fields from the finest patch level, gather from the coarsest
Definition MultiParticleContainer.H:417
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:1659
void DepositCharge(const ablastr::fields::MultiLevelScalarField &rho, amrex::Real relative_time)
Deposit charge density.
Definition MultiParticleContainer.cpp:605
void CheckQEDProductSpecies()
Definition MultiParticleContainer.cpp:1837
amrex::Real m_qed_schwinger_y_size
Definition MultiParticleContainer.H:510
amrex::Gpu::DeviceVector< amrex::ParticleReal > d_repeated_plasma_lens_strengths_E
Definition MultiParticleContainer.H:356
amrex::ParticleReal maxParticleVelocity()
Definition MultiParticleContainer.cpp:436
int nContainers() const
Definition MultiParticleContainer.H:266
std::unique_ptr< amrex::Parser > m_E_particle_from_file_parser
Definition MultiParticleContainer.H:343
void WriteHeader(std::ostream &os) const
Definition ParticleIO.cpp:242
std::string m_qed_schwinger_pos_product_name
Definition MultiParticleContainer.H:504
void Restart(const std::string &dir)
Definition ParticleIO.cpp:129
std::vector< std::string > GetSpeciesNames() const
Definition MultiParticleContainer.H:317
void ContinuousFluxInjection(amrex::Real t, amrex::Real dt) const
Definition MultiParticleContainer.cpp:925
int m_qed_schwinger_pos_product
Definition MultiParticleContainer.H:508
amrex::ParserExecutor< 1 > m_Efield_time_partparser
Definition MultiParticleContainer.H:345
amrex::MFItInfo getMFItInfo(const WarpXParticleContainer &pc_src, Args const &... pc_dsts) const noexcept
Definition MultiParticleContainer.H:423
MultiParticleContainer(amrex::AmrCore *amr_core)
Definition MultiParticleContainer.cpp:94
std::unique_ptr< amrex::Parser > m_B_particle_from_file_parser
Definition MultiParticleContainer.H:341
std::unique_ptr< amrex::Parser > m_Ez_particle_parser
Definition MultiParticleContainer.H:339
std::shared_ptr< QuantumSynchrotronEngine > m_shr_p_qs_engine
Definition MultiParticleContainer.H:445
~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:413
void SetDoBackTransformedParticles(bool do_back_transformed_particles)
Definition MultiParticleContainer.cpp:1006
PCTypes
Definition MultiParticleContainer.H:404
@ RigidInjected
Definition MultiParticleContainer.H:404
@ Physical
Definition MultiParticleContainer.H:404
@ Photon
Definition MultiParticleContainer.H:404
std::unique_ptr< amrex::Parser > m_Bx_particle_parser
Definition MultiParticleContainer.H:333
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:1673
Definition WarpXParticleContainer.H:200
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:54
amrex::Vector< ScalarField > MultiLevelScalarField
Definition MultiFabRegister.H:200
amrex::Vector< VectorField > MultiLevelVectorField
Definition MultiFabRegister.H:208
bool notInLaunchRegion() noexcept
PODVector< T, ArenaAllocator< T > > DeviceVector
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