11#ifndef WARPX_ParticleContainer_H_
12#define WARPX_ParticleContainer_H_
33#include <AMReX_Config.H>
113 std::string
const& current_fp_string,
117 bool skip_deposition=
false,
199 std::unique_ptr<amrex::MultiFab>
GetChargeDensity(
int lev,
bool local =
false);
208 void CalculateNuei (std::string
const & electron_species_name);
242 void Restart (
const std::string& dir);
252 bool sort_particles_for_deposition,
308 bool const onMainGrid =
true;
310 return static_cast<int>(std::count(
v.begin(),
v.end(), onMainGrid ));
315 bool const fromMainGrid =
true;
317 return static_cast<int>(std::count(
v.begin(),
v.end(), fromMainGrid ));
385 [[nodiscard]]
int getSpeciesID (
const std::string& product_str)
const;
434 template<
typename ...Args>
437 Args
const&... pc_dsts)
const noexcept
555 template<
typename First,
typename ...Args>
557 First
const& pc_dst, Args
const&... others)
const noexcept
561 "For particle creation processes, either all or none of the "
562 "particle species must use tiling.");
@ 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 ¤t_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_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
Definition ImplicitOptions.H:7
Definition MultiFabRegister.H:272
MFItInfo & SetDynamic(bool f) noexcept
MFItInfo & EnableTiling(const IntVect &ts=FabArrayBase::mfiter_tile_size) noexcept