WarpX
Public Member Functions | Public Attributes | Protected Types | Protected Member Functions | Protected Attributes | Static Protected Attributes | Private Member Functions | Private Attributes | List of all members
MultiParticleContainer Class Reference

#include <MultiParticleContainer.H>

Public Member Functions

 MultiParticleContainer (amrex::AmrCore *amr_core)
 
 ~MultiParticleContainer ()
 
WarpXParticleContainerGetParticleContainer (int ispecies) const
 
WarpXParticleContainerGetParticleContainerPtr (int ispecies) const
 
std::unique_ptr< WarpXParticleContainer > & GetUniqueContainer (int ispecies)
 
std::array< amrex::Real, 3 > meanParticleVelocity (int ispecies)
 
void AllocData ()
 
void InitData ()
 
void Evolve (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, const amrex::MultiFab &Ex_avg, const amrex::MultiFab &Ey_avg, const amrex::MultiFab &Ez_avg, const amrex::MultiFab &Bx_avg, const amrex::MultiFab &By_avg, const amrex::MultiFab &Bz_avg, amrex::MultiFab &jx, amrex::MultiFab &jy, amrex::MultiFab &jz, amrex::MultiFab *cjx, amrex::MultiFab *cjy, amrex::MultiFab *cjz, amrex::MultiFab *rho, amrex::MultiFab *crho, const amrex::MultiFab *cEx, const amrex::MultiFab *cEy, const amrex::MultiFab *cEz, const amrex::MultiFab *cBx, const amrex::MultiFab *cBy, const amrex::MultiFab *cBz, amrex::Real t, amrex::Real dt, DtType a_dt_type=DtType::Full)
 
void PushX (amrex::Real dt)
 
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)
 
std::unique_ptr< amrex::MultiFab > GetChargeDensity (int lev, bool local=false)
 
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)
 
void doCoulombCollisions ()
 
void doResampling (const int timestep)
 This function loops over all species and performs resampling if appropriate. More...
 
void doQEDSchwinger ()
 
void Restart (const std::string &dir)
 
void PostRestart ()
 
void ReadHeader (std::istream &is)
 
void WriteHeader (std::ostream &os) const
 
void SortParticlesByBin (amrex::IntVect bin_size)
 
void Redistribute ()
 
void RedistributeLocal (const int num_ghost)
 
void ApplyBoundaryConditions ()
 
amrex::Vector< long > NumberOfParticlesInGrid (int lev) const
 
void Increment (amrex::MultiFab &mf, int lev)
 
void SetParticleBoxArray (int lev, amrex::BoxArray &new_ba)
 
void SetParticleDistributionMap (int lev, amrex::DistributionMapping &new_dm)
 
int nSpecies () const
 
int nSpeciesBackTransformedDiagnostics () const
 
int mapSpeciesBackTransformedDiagnostics (int i) const
 
int doBackTransformedDiagnostics () const
 
int nSpeciesDepositOnMainGrid () const
 
int nSpeciesGatherFromMainGrid () const
 
void GetLabFrameData (const std::string &snapshot_name, const int i_lab, const int direction, const amrex::Real z_old, const amrex::Real z_new, const amrex::Real t_boost, const amrex::Real t_lab, const amrex::Real dt, amrex::Vector< WarpXParticleContainer::DiagnosticParticleData > &parts) const
 
void ContinuousInjection (const amrex::RealBox &injection_box) const
 
void UpdateContinuousInjectionPosition (amrex::Real dt) const
 
int doContinuousInjection () const
 
std::vector< std::string > GetSpeciesNames () const
 
PhysicalParticleContainerGetPCtmp ()
 
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) More...
 
int getSpeciesID (std::string product_str) const
 

Public Attributes

std::string m_B_ext_particle_s = "default"
 
std::string m_E_ext_particle_s = "default"
 
amrex::Vector< amrex::Real > m_B_external_particle
 
amrex::Vector< amrex::Real > m_E_external_particle
 
std::unique_ptr< ParserWrapper< 4 > > m_Bx_particle_parser
 
std::unique_ptr< ParserWrapper< 4 > > m_By_particle_parser
 
std::unique_ptr< ParserWrapper< 4 > > m_Bz_particle_parser
 
std::unique_ptr< ParserWrapper< 4 > > m_Ex_particle_parser
 
std::unique_ptr< ParserWrapper< 4 > > m_Ey_particle_parser
 
std::unique_ptr< ParserWrapper< 4 > > m_Ez_particle_parser
 

Protected Types

enum  PCTypes { PCTypes::Physical, PCTypes::RigidInjected, PCTypes::Photon }
 

Protected Member Functions

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. More...
 
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. More...
 
template<typename ... Args>
amrex::MFItInfo getMFItInfo (const WarpXParticleContainer &pc_src, Args const &... pc_dsts) const noexcept
 
void InitQED ()
 
int NSpeciesQuantumSync () const
 
int NSpeciesBreitWheeler () const
 
void InitQuantumSync ()
 
void InitBreitWheeler ()
 
void QuantumSyncGenerateTable ()
 
void BreitWheelerGenerateTable ()
 

Protected Attributes

std::vector< std::string > species_names
 
std::vector< std::string > lasers_names
 
std::vector< std::string > collision_names
 
amrex::Vector< std::unique_ptr< CollisionType > > allcollisions
 
std::vector< bool > m_deposit_on_main_grid
 instead of depositing (current, charge) on the finest patch level, deposit to the coarsest grid More...
 
std::vector< bool > m_gather_from_main_grid
 instead of gathering fields from the finest patch level, gather from the coarsest More...
 
std::vector< PCTypesspecies_types
 
ParticleBC m_boundary_conditions = ParticleBC::none
 
std::shared_ptr< BreitWheelerEnginem_shr_p_bw_engine
 
std::shared_ptr< QuantumSynchrotronEnginem_shr_p_qs_engine
 
int m_nspecies_quantum_sync = 0
 
int m_nspecies_breit_wheeler = 0
 
amrex::ParticleReal m_quantum_sync_photon_creation_energy_threshold
 
bool m_do_qed_schwinger = 0
 
std::string m_qed_schwinger_ele_product_name
 
std::string m_qed_schwinger_pos_product_name
 
int m_qed_schwinger_ele_product
 
int m_qed_schwinger_pos_product
 
amrex::Real m_qed_schwinger_y_size
 
int m_qed_schwinger_threshold_poisson_gaussian = 25
 

Static Protected Attributes

static constexpr auto m_default_quantum_sync_photon_creation_energy_threshold
 

Private Member Functions

void ReadParameters ()
 
void mapSpeciesProduct ()
 
void MFItInfoCheckTiling (const WarpXParticleContainer &) const noexcept
 
template<typename First , typename ... Args>
void MFItInfoCheckTiling (const WarpXParticleContainer &pc_src, First const &pc_dst, Args const &... others) const noexcept
 
void CheckIonizationProductSpecies ()
 
void CheckQEDProductSpecies ()
 

Private Attributes

amrex::Vector< std::unique_ptr< WarpXParticleContainer > > allcontainers
 
std::unique_ptr< PhysicalParticleContainerpc_tmp
 
int nspecies_back_transformed_diagnostics = 0
 
std::vector< int > map_species_back_transformed_diagnostics
 
int do_back_transformed_diagnostics = 0
 

Detailed Description

The class MultiParticleContainer holds multiple instances of the polymorphic class WarpXParticleContainer, stored in its member variable "allcontainers". The class WarpX typically has a single (pointer to an) instance of MultiParticleContainer.

MultiParticleContainer typically has two types of functions:

Member Enumeration Documentation

◆ PCTypes

enum MultiParticleContainer::PCTypes
strongprotected
Enumerator
Physical 
RigidInjected 
Photon 

Constructor & Destructor Documentation

◆ MultiParticleContainer()

MultiParticleContainer::MultiParticleContainer ( amrex::AmrCore *  amr_core)

◆ ~MultiParticleContainer()

MultiParticleContainer::~MultiParticleContainer ( )
inline

Member Function Documentation

◆ AllocData()

void MultiParticleContainer::AllocData ( )

◆ ApplyBoundaryConditions()

void MultiParticleContainer::ApplyBoundaryConditions ( )

Apply BC. For now, just discard particles outside the domain, regardless of the whole simulation BC.

◆ BreitWheelerGenerateTable()

void MultiParticleContainer::BreitWheelerGenerateTable ( )
protected

Called by InitBreitWheeler if a new table has to be generated.

◆ CheckIonizationProductSpecies()

void MultiParticleContainer::CheckIonizationProductSpecies ( )
private

Should be called right after mapSpeciesProduct in InitData. It checks the physical correctness of product particle species selected by the user for ionization process.

◆ CheckQEDProductSpecies()

void MultiParticleContainer::CheckQEDProductSpecies ( )
private

Should be called right after mapSpeciesProduct in InitData. It checks the physical correctness of product particle species selected by the user for QED processes.

◆ ContinuousInjection()

void MultiParticleContainer::ContinuousInjection ( const amrex::RealBox &  injection_box) const

◆ doBackTransformedDiagnostics()

int MultiParticleContainer::doBackTransformedDiagnostics ( ) const
inline

◆ doContinuousInjection()

int MultiParticleContainer::doContinuousInjection ( ) const

◆ doCoulombCollisions()

void MultiParticleContainer::doCoulombCollisions ( )

◆ doFieldIonization()

void MultiParticleContainer::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 
)

◆ doQedBreitWheeler()

void MultiParticleContainer::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 
)
protected

Performs Breit-Wheeler process for the species for which it is enabled.

◆ doQedEvents()

void MultiParticleContainer::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)

◆ doQedQuantumSync()

void MultiParticleContainer::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 
)
protected

Performs QED photon emission for the species for which it is enabled.

◆ doQEDSchwinger()

void MultiParticleContainer::doQEDSchwinger ( )

If Schwinger process is activated, this function is called at every timestep in Evolve and is used to create Schwinger electron-positron pairs. Within this function we loop over all cells to calculate the number of created physical pairs. If this number is higher than 0, we create a single particle per species in this cell, with a weight corresponding to the number of physical particles.

◆ doResampling()

void MultiParticleContainer::doResampling ( const int  timestep)

This function loops over all species and performs resampling if appropriate.

Parameters
[in]timestepthe current timestep.

◆ Evolve()

void MultiParticleContainer::Evolve ( 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,
const amrex::MultiFab &  Ex_avg,
const amrex::MultiFab &  Ey_avg,
const amrex::MultiFab &  Ez_avg,
const amrex::MultiFab &  Bx_avg,
const amrex::MultiFab &  By_avg,
const amrex::MultiFab &  Bz_avg,
amrex::MultiFab &  jx,
amrex::MultiFab &  jy,
amrex::MultiFab &  jz,
amrex::MultiFab *  cjx,
amrex::MultiFab *  cjy,
amrex::MultiFab *  cjz,
amrex::MultiFab *  rho,
amrex::MultiFab *  crho,
const amrex::MultiFab *  cEx,
const amrex::MultiFab *  cEy,
const amrex::MultiFab *  cEz,
const amrex::MultiFab *  cBx,
const amrex::MultiFab *  cBy,
const amrex::MultiFab *  cBz,
amrex::Real  t,
amrex::Real  dt,
DtType  a_dt_type = DtType::Full 
)

This evolves all the particles by one PIC time step, including current deposition, the field solve, and pushing the particles, for all the species in the MultiParticleContainer. This is the electromagnetic version.

◆ GetChargeDensity()

std::unique_ptr< MultiFab > MultiParticleContainer::GetChargeDensity ( int  lev,
bool  local = false 
)

This deposits the particle charge onto a node-centered MultiFab and returns a unique ptr to it. The charge density is accumulated over all the particles in the MultiParticleContainer

◆ GetLabFrameData()

void MultiParticleContainer::GetLabFrameData ( const std::string &  snapshot_name,
const int  i_lab,
const int  direction,
const amrex::Real  z_old,
const amrex::Real  z_new,
const amrex::Real  t_boost,
const amrex::Real  t_lab,
const amrex::Real  dt,
amrex::Vector< WarpXParticleContainer::DiagnosticParticleData > &  parts 
) const

◆ getMFItInfo()

template<typename ... Args>
amrex::MFItInfo MultiParticleContainer::getMFItInfo ( const WarpXParticleContainer pc_src,
Args const &...  pc_dsts 
) const
inlineprotectednoexcept

◆ GetParticleContainer()

WarpXParticleContainer& MultiParticleContainer::GetParticleContainer ( int  ispecies) const
inline

◆ GetParticleContainerPtr()

WarpXParticleContainer* MultiParticleContainer::GetParticleContainerPtr ( int  ispecies) const
inline

◆ GetPCtmp()

PhysicalParticleContainer& MultiParticleContainer::GetPCtmp ( )
inline

◆ getSpeciesID()

int MultiParticleContainer::getSpeciesID ( std::string  product_str) const

◆ GetSpeciesNames()

std::vector<std::string> MultiParticleContainer::GetSpeciesNames ( ) const
inline

◆ GetUniqueContainer()

std::unique_ptr<WarpXParticleContainer>& MultiParticleContainer::GetUniqueContainer ( int  ispecies)
inline

◆ Increment()

void MultiParticleContainer::Increment ( amrex::MultiFab &  mf,
int  lev 
)

◆ InitBreitWheeler()

void MultiParticleContainer::InitBreitWheeler ( )
protected

Initializes the Quantum Synchrotron engine

◆ InitData()

void MultiParticleContainer::InitData ( )

◆ InitQED()

void MultiParticleContainer::InitQED ( )
protected

Initialize QED engines and provides smart pointers to species who need QED processes

◆ InitQuantumSync()

void MultiParticleContainer::InitQuantumSync ( )
protected

Initializes the Quantum Synchrotron engine

◆ mapSpeciesBackTransformedDiagnostics()

int MultiParticleContainer::mapSpeciesBackTransformedDiagnostics ( int  i) const
inline

◆ mapSpeciesProduct()

void MultiParticleContainer::mapSpeciesProduct ( )
private

◆ meanParticleVelocity()

std::array<amrex::Real, 3> MultiParticleContainer::meanParticleVelocity ( int  ispecies)
inline

◆ MFItInfoCheckTiling() [1/2]

void MultiParticleContainer::MFItInfoCheckTiling ( const WarpXParticleContainer ) const
inlineprivatenoexcept

◆ MFItInfoCheckTiling() [2/2]

template<typename First , typename ... Args>
void MultiParticleContainer::MFItInfoCheckTiling ( const WarpXParticleContainer pc_src,
First const &  pc_dst,
Args const &...  others 
) const
inlineprivatenoexcept

◆ nSpecies()

int MultiParticleContainer::nSpecies ( ) const
inline

◆ nSpeciesBackTransformedDiagnostics()

int MultiParticleContainer::nSpeciesBackTransformedDiagnostics ( ) const
inline

◆ NSpeciesBreitWheeler()

int MultiParticleContainer::NSpeciesBreitWheeler ( ) const
inlineprotected

Returns the number of species having Breit Wheeler process enabled

◆ nSpeciesDepositOnMainGrid()

int MultiParticleContainer::nSpeciesDepositOnMainGrid ( ) const
inline

◆ nSpeciesGatherFromMainGrid()

int MultiParticleContainer::nSpeciesGatherFromMainGrid ( ) const
inline

◆ NSpeciesQuantumSync()

int MultiParticleContainer::NSpeciesQuantumSync ( ) const
inlineprotected

Returns the number of species having Quantum Synchrotron process enabled

◆ NumberOfParticlesInGrid()

Vector< long > MultiParticleContainer::NumberOfParticlesInGrid ( int  lev) const

◆ PostRestart()

void MultiParticleContainer::PostRestart ( )

◆ PushP()

void MultiParticleContainer::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 
)

This pushes the particle momenta by dt for all the species in the MultiParticleContainer. It is used to desynchronize the particles after initializaton or when restarting from a checkpoint. It is also used to synchronize particles at the the end of the run. This is the electromagnetic version.

◆ PushX()

void MultiParticleContainer::PushX ( amrex::Real  dt)

This pushes the particle positions by one half time step for all the species in the MultiParticleContainer. It is used to desynchronize the particles after initializaton or when restarting from a checkpoint.

◆ QuantumSyncGenerateTable()

void MultiParticleContainer::QuantumSyncGenerateTable ( )
protected

Called by InitQuantumSync if a new table has to be generated.

◆ ReadHeader()

void MultiParticleContainer::ReadHeader ( std::istream &  is)

◆ ReadParameters()

void MultiParticleContainer::ReadParameters ( )
private

◆ Redistribute()

void MultiParticleContainer::Redistribute ( )

◆ RedistributeLocal()

void MultiParticleContainer::RedistributeLocal ( const int  num_ghost)

◆ Restart()

void MultiParticleContainer::Restart ( const std::string &  dir)

◆ SetParticleBoxArray()

void MultiParticleContainer::SetParticleBoxArray ( int  lev,
amrex::BoxArray &  new_ba 
)

◆ SetParticleDistributionMap()

void MultiParticleContainer::SetParticleDistributionMap ( int  lev,
amrex::DistributionMapping &  new_dm 
)

◆ SortParticlesByBin()

void MultiParticleContainer::SortParticlesByBin ( amrex::IntVect  bin_size)

◆ UpdateContinuousInjectionPosition()

void MultiParticleContainer::UpdateContinuousInjectionPosition ( amrex::Real  dt) const

◆ WriteHeader()

void MultiParticleContainer::WriteHeader ( std::ostream &  os) const

Member Data Documentation

◆ allcollisions

amrex::Vector<std::unique_ptr<CollisionType> > MultiParticleContainer::allcollisions
protected

◆ allcontainers

amrex::Vector<std::unique_ptr<WarpXParticleContainer> > MultiParticleContainer::allcontainers
private

◆ collision_names

std::vector<std::string> MultiParticleContainer::collision_names
protected

◆ do_back_transformed_diagnostics

int MultiParticleContainer::do_back_transformed_diagnostics = 0
private

◆ lasers_names

std::vector<std::string> MultiParticleContainer::lasers_names
protected

◆ m_B_ext_particle_s

std::string MultiParticleContainer::m_B_ext_particle_s = "default"

◆ m_B_external_particle

amrex::Vector<amrex::Real> MultiParticleContainer::m_B_external_particle

◆ m_boundary_conditions

ParticleBC MultiParticleContainer::m_boundary_conditions = ParticleBC::none
protected

Whether to absorb particles exiting the domain

◆ m_Bx_particle_parser

std::unique_ptr<ParserWrapper<4> > MultiParticleContainer::m_Bx_particle_parser

◆ m_By_particle_parser

std::unique_ptr<ParserWrapper<4> > MultiParticleContainer::m_By_particle_parser

◆ m_Bz_particle_parser

std::unique_ptr<ParserWrapper<4> > MultiParticleContainer::m_Bz_particle_parser

◆ m_default_quantum_sync_photon_creation_energy_threshold

constexpr auto MultiParticleContainer::m_default_quantum_sync_photon_creation_energy_threshold
staticprotected
Initial value:
=
static_cast<amrex::ParticleReal>(
2.0 * PhysConst::m_e * PhysConst::c * PhysConst::c )

Default value of the energy threshold for photon creation in Quantum Synchrotron process.

◆ m_deposit_on_main_grid

std::vector<bool> MultiParticleContainer::m_deposit_on_main_grid
protected

instead of depositing (current, charge) on the finest patch level, deposit to the coarsest grid

◆ m_do_qed_schwinger

bool MultiParticleContainer::m_do_qed_schwinger = 0
protected

Whether or not to activate Schwinger process

◆ m_E_ext_particle_s

std::string MultiParticleContainer::m_E_ext_particle_s = "default"

◆ m_E_external_particle

amrex::Vector<amrex::Real> MultiParticleContainer::m_E_external_particle

◆ m_Ex_particle_parser

std::unique_ptr<ParserWrapper<4> > MultiParticleContainer::m_Ex_particle_parser

◆ m_Ey_particle_parser

std::unique_ptr<ParserWrapper<4> > MultiParticleContainer::m_Ey_particle_parser

◆ m_Ez_particle_parser

std::unique_ptr<ParserWrapper<4> > MultiParticleContainer::m_Ez_particle_parser

◆ m_gather_from_main_grid

std::vector<bool> MultiParticleContainer::m_gather_from_main_grid
protected

instead of gathering fields from the finest patch level, gather from the coarsest

◆ m_nspecies_breit_wheeler

int MultiParticleContainer::m_nspecies_breit_wheeler = 0
protected

◆ m_nspecies_quantum_sync

int MultiParticleContainer::m_nspecies_quantum_sync = 0
protected

◆ m_qed_schwinger_ele_product

int MultiParticleContainer::m_qed_schwinger_ele_product
protected

Index for Schwinger electron product species in allcontainers

◆ m_qed_schwinger_ele_product_name

std::string MultiParticleContainer::m_qed_schwinger_ele_product_name
protected

Name of Schwinger electron product species

◆ m_qed_schwinger_pos_product

int MultiParticleContainer::m_qed_schwinger_pos_product
protected

Index for Schwinger positron product species in allcontainers

◆ m_qed_schwinger_pos_product_name

std::string MultiParticleContainer::m_qed_schwinger_pos_product_name
protected

Name of Schwinger positron product species

◆ m_qed_schwinger_threshold_poisson_gaussian

int MultiParticleContainer::m_qed_schwinger_threshold_poisson_gaussian = 25
protected

If the number of physical Schwinger pairs created within a cell is higher than this threshold we use a Gaussian distribution rather than a Poisson distribution for the pair production rate calculations

◆ m_qed_schwinger_y_size

amrex::Real MultiParticleContainer::m_qed_schwinger_y_size
protected

Transverse size used in 2D Schwinger pair production rate calculations

◆ m_quantum_sync_photon_creation_energy_threshold

amrex::ParticleReal MultiParticleContainer::m_quantum_sync_photon_creation_energy_threshold
protected
Initial value:

Energy threshold for photon creation in Quantum Synchrotron process.

◆ m_shr_p_bw_engine

std::shared_ptr<BreitWheelerEngine> MultiParticleContainer::m_shr_p_bw_engine
protected

◆ m_shr_p_qs_engine

std::shared_ptr<QuantumSynchrotronEngine> MultiParticleContainer::m_shr_p_qs_engine
protected

◆ map_species_back_transformed_diagnostics

std::vector<int> MultiParticleContainer::map_species_back_transformed_diagnostics
private

◆ nspecies_back_transformed_diagnostics

int MultiParticleContainer::nspecies_back_transformed_diagnostics = 0
private

◆ pc_tmp

std::unique_ptr<PhysicalParticleContainer> MultiParticleContainer::pc_tmp
private

◆ species_names

std::vector<std::string> MultiParticleContainer::species_names
protected

◆ species_types

std::vector<PCTypes> MultiParticleContainer::species_types
protected

The documentation for this class was generated from the following files: