ImpactX
Loading...
Searching...
No Matches
impactx::particles::spacecharge Namespace Reference

Functions

void ForceFromSelfFields (std::unordered_map< int, std::unordered_map< std::string, amrex::MultiFab > > &space_charge_field, std::unordered_map< int, amrex::MultiFab > const &phi, amrex::Vector< amrex::Geometry > const &geom)
 
void GatherAndPush (ImpactXParticleContainer &pc, std::unordered_map< int, std::unordered_map< std::string, amrex::MultiFab > > const &space_charge_field, const amrex::Vector< amrex::Geometry > &geom, amrex::ParticleReal const slice_ds)
 
AMREX_GPU_DEVICE void efldgauss (int nint, amrex::ParticleReal const xin, amrex::ParticleReal const yin, amrex::ParticleReal const zin, amrex::ParticleReal const sigx, amrex::ParticleReal const sigy, amrex::ParticleReal const sigz, amrex::ParticleReal const gamma, amrex::ParticleReal &pintex, amrex::ParticleReal &pintey, amrex::ParticleReal &pintez)
 
void Gauss3dPush (ImpactXParticleContainer &pc, amrex::ParticleReal const slice_ds)
 
void HandleSpacecharge (std::unique_ptr< initialization::AmrCoreData > &amr_data, std::function< void()> ResizeMesh, amrex::Real slice_ds)
 
void PoissonSolve (ImpactXParticleContainer const &pc, std::unordered_map< int, amrex::MultiFab > &rho, std::unordered_map< int, amrex::MultiFab > &phi, amrex::Vector< amrex::IntVect > rel_ref_ratio)
 

Function Documentation

◆ efldgauss()

AMREX_GPU_DEVICE void impactx::particles::spacecharge::efldgauss ( int nint,
amrex::ParticleReal const xin,
amrex::ParticleReal const yin,
amrex::ParticleReal const zin,
amrex::ParticleReal const sigx,
amrex::ParticleReal const sigy,
amrex::ParticleReal const sigz,
amrex::ParticleReal const gamma,
amrex::ParticleReal & pintex,
amrex::ParticleReal & pintey,
amrex::ParticleReal & pintez )

Compute space-charge fields from a 3D Gaussian distribution.

Compute integrals Eqs 45-47 used in the space-charge fields from a 3D Gaussian distribution. "A two-and-a-half dimensional symplectic space-charge solver", LBNL Report Number: LBNL-2001674, 2025. Input particle locations (x,y,z) and RMS sizes (sigx,sigy,sigz) and return the integrals for SC fields.

Todo
This function can be vectorized with AMREX_SIMD for better CPU performance.

◆ ForceFromSelfFields()

void impactx::particles::spacecharge::ForceFromSelfFields ( std::unordered_map< int, std::unordered_map< std::string, amrex::MultiFab > > & space_charge_field,
std::unordered_map< int, amrex::MultiFab > const & phi,
const amrex::Vector< amrex::Geometry > & geom )

Calculate the space charge force field from the electric potential

This resets the values in scf_<component> to zero and then calculates the space charge force field.

Parameters
[in,out]space_charge_fieldspace charge force component in x,y,z per level
[in]phiscalar potential per level
[in]geomgeometry object

◆ GatherAndPush()

void impactx::particles::spacecharge::GatherAndPush ( ImpactXParticleContainer & pc,
std::unordered_map< int, std::unordered_map< std::string, amrex::MultiFab > > const & space_charge_field,
const amrex::Vector< amrex::Geometry > & geom,
amrex::ParticleReal slice_ds )

Gather force fields and push particles in x,y,z

This gathers the space charge field with respect to particle position and shape. The momentum of all particles is then pushed using a common time step given by the reference particle speed and ds slice. The position push is done in the lattice elements and not here.

Parameters
[in,out]pccontainer of the particles that deposited rho
[in]space_charge_fieldspace charge force component in x,y,z per level
[in]geomgeometry object
[in]slice_dssegment length in meters

◆ Gauss3dPush()

void impactx::particles::spacecharge::Gauss3dPush ( ImpactXParticleContainer & pc,
amrex::ParticleReal slice_ds )

Apply a Space-Charge Kick according to a 3D Gaussian Distribution

This calculates the space charge fields from a 3D Gaussian distribution with respect to particle position from the following paper: J. Qiang et al., "Two-and-a-half dimensional symplectic space-charge solver", LBNL Report Number: LBNL-2001674 (2025). The momentum of all particles is then pushed using a common time step given by the reference particle speed and ds slice. The position push is done in the lattice elements and not here.

Parameters
[in,out]pccontainer of the particles that deposited rho
[in]slice_dssegment length in meters

◆ HandleSpacecharge()

void impactx::particles::spacecharge::HandleSpacecharge ( std::unique_ptr< initialization::AmrCoreData > & amr_data,
std::function< void()> ResizeMesh,
amrex::Real slice_ds )

Function to handle space charge including charge deposition, Poisson solve, field calculation, interpolation, and particle push.

Parameters
[in,out]amr_dataAMR data like particle container and fields
[in]ResizeMeshfunction to call to resize the mesh
[in]slice_dsslice spacing along s

◆ PoissonSolve()

void impactx::particles::spacecharge::PoissonSolve ( ImpactXParticleContainer const & pc,
std::unordered_map< int, amrex::MultiFab > & rho,
std::unordered_map< int, amrex::MultiFab > & phi,
amrex::Vector< amrex::IntVect > rel_ref_ratio )

Calculate the electric potential from charge density

This resets the values in phi to zero and then calculates the space charge potential phi.

Parameters
[in]pccontainer of the particles that deposited rho
[in]rhocharge per level
[in,out]phiscalar potential per level
[in]rel_ref_ratiomesh refinement ratio between levels