|
ImpactX
|
| amrex::Gpu::DeviceVector< amrex::Real > impactx::particles::spacecharge::Deposit1D | ( | ImpactXParticleContainer & | pc, |
| amrex::Real | bin_min, | ||
| amrex::Real | bin_max, | ||
| int | num_bins ) |
Deposit charge using longitudinal binning for the calculation of space charge
This function performs particle binning for particles in the beam along the longitudinal direction, for the purposes of computing space charge in so-called 2.5D models.
| [in] | pc | container of the particles that deposited rho |
| [in] | bin_min | location of lower endpoint for binning |
| [in] | bin_max | location of upper endpoint for binning |
| [in] | num_bins | number of longitudinal bins |
| 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.
| 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.
| [in,out] | space_charge_field | space charge force component in x,y,z per level |
| [in] | phi | scalar potential per level |
| [in] | geom | geometry object |
| void impactx::particles::spacecharge::GatherAndPush | ( | ImpactXParticleContainer & | pc, |
| std::unordered_map< int, std::unordered_map< std::string, amrex::MultiFab > > const & | space_charge_field, | ||
| std::unordered_map< int, amrex::MultiFab > const & | space_charge_potential, | ||
| 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.
| [in,out] | pc | container of the particles that deposited rho |
| [in] | space_charge_field | space charge force component in x,y,z per level |
| [in] | space_charge_potential | space charge potential in x,y,z per level |
| [in] | geom | geometry object |
| [in] | slice_ds | segment length in meters |
| void impactx::particles::spacecharge::Gauss2p5dPush | ( | ImpactXParticleContainer & | pc, |
| amrex::ParticleReal | slice_ds ) |
Calculate the 2.5D Space Charge fields for a transverse Gaussian Distribution
This calculates the space charge fields from a transverse Gaussian distribution with respect to particle position from the following paper. "Two-and-a-half dimensional symplectic space-charge solver" 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.
| [in,out] | pc | container of the particles that deposited rho |
| [in] | slice_ds | segment length in meters |
| 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.
| [in,out] | pc | container of the particles that deposited rho |
| [in] | slice_ds | segment length in meters |
| 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.
| [in,out] | amr_data | AMR data like particle container and fields |
| [in] | ResizeMesh | function to call to resize the mesh |
| [in] | slice_ds | slice spacing along s |
| 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.
| [in] | pc | container of the particles that deposited rho |
| [in] | rho | charge per level |
| [in,out] | phi | scalar potential per level |
| [in] | rel_ref_ratio | mesh refinement ratio between levels |
| AMREX_GPU_DEVICE void impactx::particles::spacecharge::potInt | ( | amrex::ParticleReal | delta, |
| int | nint, | ||
| amrex::ParticleReal | xin, | ||
| amrex::ParticleReal | yin, | ||
| amrex::ParticleReal | sigx, | ||
| amrex::ParticleReal | sigy, | ||
| amrex::ParticleReal & | pintex, | ||
| amrex::ParticleReal & | pintey, | ||
| amrex::ParticleReal & | pintez ) |
Calculate the Integrals for Space Charge Fields
Compute integrals Eqs 25, 31,32 used in the space-charge fields from a 2D transverse Gaussian distribution (including 2A in Eq. 25). Input particle locations (x,y) and RMS sizes (sigx,sigy) and return the integrals for SC fields.