WarpX
Public Member Functions | Private Attributes | List of all members
RigidInjectedParticleContainer Class Reference

#include <RigidInjectedParticleContainer.H>

Inheritance diagram for RigidInjectedParticleContainer:
PhysicalParticleContainer WarpXParticleContainer

Public Member Functions

 RigidInjectedParticleContainer (amrex::AmrCore *amr_core, int ispecies, const std::string &name)
 
virtual ~RigidInjectedParticleContainer ()
 
virtual void InitData () override
 
virtual void RemapParticles ()
 
virtual void BoostandRemapParticles ()
 
virtual 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) override
 Evolve is the central function PhysicalParticleContainer that advances plasma particles for a time dt (typically one timestep). More...
 
virtual void PushPX (WarpXParIter &pti, amrex::FArrayBox const *exfab, amrex::FArrayBox const *eyfab, amrex::FArrayBox const *ezfab, amrex::FArrayBox const *bxfab, amrex::FArrayBox const *byfab, amrex::FArrayBox const *bzfab, const int ngE, const int, const long offset, const long np_to_push, int lev, int gather_lev, amrex::Real dt, ScaleFields scaleFields, DtType a_dt_type=DtType::Full) override
 
virtual 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) override
 
virtual void ReadHeader (std::istream &is) override
 
virtual void WriteHeader (std::ostream &os) const override
 
- Public Member Functions inherited from PhysicalParticleContainer
 PhysicalParticleContainer (amrex::AmrCore *amr_core, int ispecies, const std::string &name)
 
 PhysicalParticleContainer (amrex::AmrCore *amr_core)
 
void BackwardCompatibility ()
 
virtual ~PhysicalParticleContainer ()
 
void InitIonizationModule ()
 
void PartitionParticlesInBuffers (long &nfine_current, long &nfine_gather, long const np, WarpXParIter &pti, int const lev, amrex::iMultiFab const *current_masks, amrex::iMultiFab const *gather_masks, RealVector &uxp, RealVector &uyp, RealVector &uzp, RealVector &wp)
 
virtual void PostRestart () final
 
void SplitParticles (int lev)
 
IonizationFilterFunc getIonizationFunc (const WarpXParIter &pti, int lev, int ngE, const amrex::FArrayBox &Ex, const amrex::FArrayBox &Ey, const amrex::FArrayBox &Ez, const amrex::FArrayBox &Bx, const amrex::FArrayBox &By, const amrex::FArrayBox &Bz)
 
virtual void AddParticles (int lev)
 
void AddPlasma (int lev, amrex::RealBox part_realbox=amrex::RealBox())
 
void MapParticletoBoostedFrame (amrex::Real &x, amrex::Real &y, amrex::Real &z, std::array< amrex::Real, 3 > &u)
 
void AddGaussianBeam (const amrex::Real x_m, const amrex::Real y_m, const amrex::Real z_m, const amrex::Real x_rms, const amrex::Real y_rms, const amrex::Real z_rms, const amrex::Real x_cut, const amrex::Real y_cut, const amrex::Real z_cut, const amrex::Real q_tot, long npart, const int do_symmetrize)
 
void AddPlasmaFromFile (amrex::ParticleReal q_tot, amrex::ParticleReal z_shift)
 
void CheckAndAddParticle (amrex::Real x, amrex::Real y, amrex::Real z, std::array< amrex::Real, 3 > u, amrex::Real weight, amrex::Gpu::HostVector< amrex::ParticleReal > &particle_x, amrex::Gpu::HostVector< amrex::ParticleReal > &particle_y, amrex::Gpu::HostVector< amrex::ParticleReal > &particle_z, amrex::Gpu::HostVector< amrex::ParticleReal > &particle_ux, amrex::Gpu::HostVector< amrex::ParticleReal > &particle_uy, amrex::Gpu::HostVector< amrex::ParticleReal > &particle_uz, amrex::Gpu::HostVector< amrex::ParticleReal > &particle_w)
 
virtual void GetParticleSlice (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, DiagnosticParticles &diagnostic_particles) final
 
virtual void ConvertUnits (ConvertDirection convert_dir) override
 
void applyNCIFilter (int lev, const amrex::Box &box, amrex::Elixir &exeli, amrex::Elixir &eyeli, amrex::Elixir &ezeli, amrex::Elixir &bxeli, amrex::Elixir &byeli, amrex::Elixir &bzeli, amrex::FArrayBox &filtered_Ex, amrex::FArrayBox &filtered_Ey, amrex::FArrayBox &filtered_Ez, amrex::FArrayBox &filtered_Bx, amrex::FArrayBox &filtered_By, amrex::FArrayBox &filtered_Bz, const amrex::FArrayBox &Ex, const amrex::FArrayBox &Ey, const amrex::FArrayBox &Ez, const amrex::FArrayBox &Bx, const amrex::FArrayBox &By, const amrex::FArrayBox &Bz, amrex::FArrayBox const *&exfab, amrex::FArrayBox const *&eyfab, amrex::FArrayBox const *&ezfab, amrex::FArrayBox const *&bxfab, amrex::FArrayBox const *&byfab, amrex::FArrayBox const *&bzfab)
 Apply NCI Godfrey filter to all components of E and B before gather. More...
 
void resample (const int timestep) override final
 This function determines if resampling should be done for the current species, and if so, performs the resampling. More...
 
bool has_quantum_sync () const override
 
bool has_breit_wheeler () const override
 
void set_breit_wheeler_engine_ptr (std::shared_ptr< BreitWheelerEngine > ptr) override
 
void set_quantum_sync_engine_ptr (std::shared_ptr< QuantumSynchrotronEngine > ptr) override
 
PhotonEmissionFilterFunc getPhotonEmissionFilterFunc ()
 
PairGenerationFilterFunc getPairGenerationFilterFunc ()
 
- Public Member Functions inherited from WarpXParticleContainer
 WarpXParticleContainer (amrex::AmrCore *amr_core, int ispecies)
 
virtual ~WarpXParticleContainer ()
 
void AllocData ()
 
void PushX (amrex::Real dt)
 
void PushX (int lev, amrex::Real dt)
 
void DepositCharge (amrex::Vector< std::unique_ptr< amrex::MultiFab > > &rho, bool local=false, bool reset=false, bool do_rz_volume_scaling=false)
 
std::unique_ptr< amrex::MultiFab > GetChargeDensity (int lev, bool local=false)
 
virtual void DepositCharge (WarpXParIter &pti, RealVector &wp, const int *const ion_lev, amrex::MultiFab *rho, int icomp, const long offset, const long np_to_depose, int thread_num, int lev, int depos_lev)
 
virtual void DepositCurrent (WarpXParIter &pti, RealVector &wp, RealVector &uxp, RealVector &uyp, RealVector &uzp, const int *const ion_lev, amrex::MultiFab *jx, amrex::MultiFab *jy, amrex::MultiFab *jz, const long offset, const long np_to_depose, int thread_num, int lev, int depos_lev, amrex::Real dt)
 
virtual void UpdateContinuousInjectionPosition (amrex::Real)
 
amrex::Real sumParticleCharge (bool local=false)
 
std::array< amrex::Real, 3 > meanParticleVelocity (bool local=false)
 
amrex::Real maxParticleVelocity (bool local=false)
 
void AddNParticles (int lev, int n, const amrex::ParticleReal *x, const amrex::ParticleReal *y, const amrex::ParticleReal *z, const amrex::ParticleReal *vx, const amrex::ParticleReal *vy, const amrex::ParticleReal *vz, int nattr, const amrex::ParticleReal *attr, int uniqueparticles, amrex::Long id=-1)
 
void ApplyBoundaryConditions (ParticleBC boundary_conditions)
 Apply particle BC. More...
 
void AddRealComp (const std::string &name, bool comm=true)
 
void AddIntComp (const std::string &name, bool comm=true)
 
int doBackTransformedDiagnostics () const
 
std::map< std::string, int > getParticleComps () const noexcept
 
std::map< std::string, int > getParticleiComps () const noexcept
 
std::map< std::string, int > getParticleRuntimeComps () const noexcept
 
std::map< std::string, int > getParticleRuntimeiComps () const noexcept
 
amrex::ParticleReal getCharge () const
 
amrex::ParticleReal getMass () const
 
int DoFieldIonization () const
 
int DoQED () const
 
template<PhysicalSpecies PhysSpec>
bool AmIA () const noexcept
 
amrex::Array< amrex::Real, 3 > get_v_galilean ()
 

Private Attributes

amrex::Real zinject_plane = 0.
 
bool projected = true
 
bool focused = false
 
bool rigid_advance = true
 
amrex::Real vzbeam_ave_boosted
 
amrex::Vector< int > done_injecting
 
amrex::Vector< amrex::Real > zinject_plane_levels
 
amrex::Real zinject_plane_lev
 
amrex::Real zinject_plane_lev_previous
 
bool done_injecting_lev
 

Additional Inherited Members

- Public Types inherited from PhysicalParticleContainer
enum  PhysicalParticleType { electron, positron, photon, other }
 
- Public Types inherited from WarpXParticleContainer
using DiagnosticParticleData = amrex::StructOfArrays< DiagIdx::nattribs, 0 >
 
using DiagnosticParticles = amrex::Vector< std::map< std::pair< int, int >, DiagnosticParticleData > >
 
using PairIndex = std::pair< int, int >
 
using TmpParticleTile = std::array< amrex::Gpu::DeviceVector< amrex::ParticleReal >, TmpIdx::nattribs >
 
using TmpParticles = amrex::Vector< std::map< PairIndex, TmpParticleTile > >
 
- Static Public Member Functions inherited from WarpXParticleContainer
static void ReadParameters ()
 
static void BackwardCompatibility ()
 
- Public Attributes inherited from WarpXParticleContainer
friend MultiParticleContainer
 
bool do_splitting = false
 
bool initialize_self_fields = false
 
amrex::Real self_fields_required_precision = 1.e-11
 
int split_type = 0
 
bool m_do_random_filter = false
 
bool m_do_uniform_filter = false
 
bool m_do_parser_filter = false
 
amrex::Real m_random_fraction = 1.0
 
int m_uniform_stride = 1
 
std::unique_ptr< ParserWrapper< 7 > > m_particle_filter_parser
 
- Protected Member Functions inherited from PhysicalParticleContainer
void ContinuousInjection (const amrex::RealBox &injection_box) override
 
- Protected Member Functions inherited from WarpXParticleContainer
void defineAllParticleTiles () noexcept
 
- Protected Attributes inherited from PhysicalParticleContainer
std::string species_name
 
std::unique_ptr< PlasmaInjectorplasma_injector
 
bool boost_adjust_transverse_positions = false
 
bool do_backward_propagation = false
 
Resampling m_resampler
 
bool do_classical_radiation_reaction = false
 
bool m_do_qed_quantum_sync = false
 
bool m_do_qed_breit_wheeler = false
 
std::shared_ptr< QuantumSynchrotronEnginem_shr_p_qs_engine
 
std::shared_ptr< BreitWheelerEnginem_shr_p_bw_engine
 
- Protected Attributes inherited from WarpXParticleContainer
amrex::Array< amrex::Real, 3 > m_v_galilean = {{0}}
 
std::map< std::string, int > particle_comps
 
std::map< std::string, int > particle_icomps
 
std::map< std::string, int > particle_runtime_comps
 
std::map< std::string, int > particle_runtime_icomps
 
int species_id
 
amrex::Real charge
 
amrex::Real mass
 
PhysicalSpecies physical_species
 
bool m_deposit_on_main_grid = false
 instead of depositing (current, charge) on the finest patch level, deposit to the coarsest grid More...
 
bool m_gather_from_main_grid = false
 instead of gathering fields from the finest patch level, gather from the coarsest More...
 
int do_not_push = 0
 
int do_not_deposit = 0
 
int do_not_gather = 0
 
int do_continuous_injection = 0
 
int do_field_ionization = 0
 
int ionization_product
 
std::string ionization_product_name
 
int ion_atomic_number
 
int ionization_initial_level = 0
 
amrex::Gpu::DeviceVector< amrex::Real > ionization_energies
 
amrex::Gpu::DeviceVector< amrex::Real > adk_power
 
amrex::Gpu::DeviceVector< amrex::Real > adk_prefactor
 
amrex::Gpu::DeviceVector< amrex::Real > adk_exp_prefactor
 
std::string physical_element
 
int do_resampling = 0
 
int do_back_transformed_diagnostics = 1
 
bool m_do_qed = false
 
int m_qed_breit_wheeler_ele_product
 
std::string m_qed_breit_wheeler_ele_product_name
 
int m_qed_breit_wheeler_pos_product
 
std::string m_qed_breit_wheeler_pos_product_name
 
int m_qed_quantum_sync_phot_product
 
std::string m_qed_quantum_sync_phot_product_name
 
amrex::Vector< amrex::FArrayBox > local_rho
 
amrex::Vector< amrex::FArrayBox > local_jx
 
amrex::Vector< amrex::FArrayBox > local_jy
 
amrex::Vector< amrex::FArrayBox > local_jz
 
TmpParticles tmp_particle_data
 

Detailed Description

When injecting a particle beam (typically for a plasma wakefield acceleration simulation), say propagating in the z direction, it can be necessary to make particles propagate in a straight line up to a given location z=z0. This is of particular importance when running in a boosted frame, where the beam may evolve due to its space charge fields before entering the plasma, causing the actual injected beam, and hence the whole simulation result, to depend on the Lorentz factor of the boost.

This feature is implemented in RigidInjectedParticleContainer: At each iteration, for each particle, if z<z0 the particle moves in a straight line, and if z>z0 the particle evolves as a regular PhysicalParticleContainer.

Note: This option is also useful to build self-consistent space charge fields for the particle beam.

RigidInjectedParticleContainer derives from PhysicalParticleContainer.

Constructor & Destructor Documentation

◆ RigidInjectedParticleContainer()

RigidInjectedParticleContainer::RigidInjectedParticleContainer ( amrex::AmrCore *  amr_core,
int  ispecies,
const std::string &  name 
)

◆ ~RigidInjectedParticleContainer()

virtual RigidInjectedParticleContainer::~RigidInjectedParticleContainer ( )
inlinevirtual

Member Function Documentation

◆ BoostandRemapParticles()

void RigidInjectedParticleContainer::BoostandRemapParticles ( )
virtual

◆ Evolve()

void RigidInjectedParticleContainer::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 
)
overridevirtual

Evolve is the central function PhysicalParticleContainer that advances plasma particles for a time dt (typically one timestep).

Parameters
levlevel on which particles are living
ExMultiFab from which field Ex is gathered
EyMultiFab from which field Ey is gathered
EzMultiFab from which field Ez is gathered
BxMultiFab from which field Bx is gathered
ByMultiFab from which field By is gathered
BzMultiFab from which field Bz is gathered
jxMultiFab to which the particles' current jx is deposited
jyMultiFab to which the particles' current jy is deposited
jzMultiFab to which the particles' current jz is deposited
cjxSame as jx (coarser, from lev-1), when using deposition buffers
cjySame as jy (coarser, from lev-1), when using deposition buffers
cjzSame as jz (coarser, from lev-1), when using deposition buffers
rhoMultiFab to which the particles' charge is deposited
crhoSame as rho (coarser, from lev-1), when using deposition buffers
cExSame as Ex (coarser, from lev-1), when using gather buffers
cEySame as Ey (coarser, from lev-1), when using gather buffers
cEzSame as Ez (coarser, from lev-1), when using gather buffers
cBxSame as Bx (coarser, from lev-1), when using gather buffers
cBySame as By (coarser, from lev-1), when using gather buffers
cBzSame as Bz (coarser, from lev-1), when using gather buffers
tcurrent physical time
dttime step by which particles are advanced
a_dt_typetype of time step (used for sub-cycling)

Evolve iterates over particle iterator (each box) and performs filtering, field gather, particle push and current deposition for all particles in the box.

Reimplemented from PhysicalParticleContainer.

◆ InitData()

void RigidInjectedParticleContainer::InitData ( )
overridevirtual

Reimplemented from PhysicalParticleContainer.

◆ PushP()

void RigidInjectedParticleContainer::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 
)
overridevirtual

This pushes the particle momenta by dt.

Reimplemented from PhysicalParticleContainer.

◆ PushPX()

void RigidInjectedParticleContainer::PushPX ( WarpXParIter pti,
amrex::FArrayBox const *  exfab,
amrex::FArrayBox const *  eyfab,
amrex::FArrayBox const *  ezfab,
amrex::FArrayBox const *  bxfab,
amrex::FArrayBox const *  byfab,
amrex::FArrayBox const *  bzfab,
const int  ngE,
const int  e_is_nodal,
const long  offset,
const long  np_to_push,
int  lev,
int  gather_lev,
amrex::Real  dt,
ScaleFields  scaleFields,
DtType  a_dt_type = DtType::Full 
)
overridevirtual

Reimplemented from PhysicalParticleContainer.

◆ ReadHeader()

void RigidInjectedParticleContainer::ReadHeader ( std::istream &  is)
overridevirtual

Reimplemented from WarpXParticleContainer.

◆ RemapParticles()

void RigidInjectedParticleContainer::RemapParticles ( )
virtual

◆ WriteHeader()

void RigidInjectedParticleContainer::WriteHeader ( std::ostream &  os) const
overridevirtual

Reimplemented from WarpXParticleContainer.

Member Data Documentation

◆ done_injecting

amrex::Vector<int> RigidInjectedParticleContainer::done_injecting
private

◆ done_injecting_lev

bool RigidInjectedParticleContainer::done_injecting_lev
private

◆ focused

bool RigidInjectedParticleContainer::focused = false
private

◆ projected

bool RigidInjectedParticleContainer::projected = true
private

◆ rigid_advance

bool RigidInjectedParticleContainer::rigid_advance = true
private

◆ vzbeam_ave_boosted

amrex::Real RigidInjectedParticleContainer::vzbeam_ave_boosted
private

◆ zinject_plane

amrex::Real RigidInjectedParticleContainer::zinject_plane = 0.
private

◆ zinject_plane_lev

amrex::Real RigidInjectedParticleContainer::zinject_plane_lev
private

◆ zinject_plane_lev_previous

amrex::Real RigidInjectedParticleContainer::zinject_plane_lev_previous
private

◆ zinject_plane_levels

amrex::Vector<amrex::Real> RigidInjectedParticleContainer::zinject_plane_levels
private

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