|
ImpactX
|
Functions | |
| void | DepositCharge1D (impactx::ImpactXParticleContainer &myspc, amrex::Gpu::DeviceVector< amrex::Real > &charge_distribution, amrex::Real bin_min, amrex::Real bin_size, bool is_unity_particle_weight) |
| void | DerivativeCharge1D (amrex::Gpu::DeviceVector< amrex::Real > const &charge_distribution, amrex::Gpu::DeviceVector< amrex::Real > &slopes, amrex::Real bin_size, bool GetNumberDensity) |
| void | MeanTransversePosition (impactx::ImpactXParticleContainer &myspc, amrex::Gpu::DeviceVector< amrex::Real > &mean_x, amrex::Gpu::DeviceVector< amrex::Real > &mean_y, amrex::Real bin_min, amrex::Real bin_size, bool is_unity_particle_weight) |
| template<typename VariantType> | |
| std::tuple< bool, amrex::Real > | CSRBendElement (VariantType const &element_variant, RefPart const &refpart) |
| template<typename T_Element> | |
| void | HandleISR (impactx::ImpactXParticleContainer &particle_container, T_Element const &element_variant, amrex::Real slice_ds) |
| template<typename T_Element> | |
| void | HandleWakefield (impactx::ImpactXParticleContainer &particle_container, T_Element const &element_variant, amrex::Real slice_ds, bool print_wakefield=false) |
| template<typename VariantType> | |
| std::tuple< bool, amrex::Real > | ISRBendElement (VariantType const &element_variant, RefPart const &refpart) |
| void | ISRPush (ImpactXParticleContainer &pc, amrex::ParticleReal slice_ds, amrex::ParticleReal rc, int isr_order) |
| amrex::Real | alpha (amrex::Real s) |
| amrex::Real | w_t_rf (amrex::Real s, amrex::Real a, amrex::Real g, amrex::Real L) |
| amrex::Real | w_l_rf (amrex::Real s, amrex::Real a, amrex::Real g, amrex::Real L) |
| amrex::Gpu::DeviceVector< amrex::Real > | convolve_fft (amrex::Gpu::DeviceVector< amrex::Real > const &beam_profile_slope, amrex::Gpu::DeviceVector< amrex::Real > const &wake_func, amrex::Real delta_t) |
| AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real | unit_step (amrex::Real s) |
| AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real | w_l_csr (amrex::Real s, amrex::Real R, amrex::Real const bin_size) |
| void | WakePush (ImpactXParticleContainer &pc, amrex::Gpu::DeviceVector< amrex::Real > const &convolved_wakefield, amrex::ParticleReal slice_ds, amrex::Real bin_size, amrex::Real bin_min) |
Variables | |
| constexpr amrex::Real | Z0 = 377 |
| Free space impedance [Ohm]. | |
| constexpr amrex::Real | alpha_1 = 0.4648 |
| amrex::Real impactx::particles::wakefields::alpha | ( | amrex::Real | s | ) |
Alpha Function
Eq. (21) in: K. L.F. Bane, "Short-Range Dipole Wakefields in Accelerating Structures for the NLC," SLAC-PUB-9663, 2003
| s | value along s [m] |
| amrex::Gpu::DeviceVector< amrex::Real > impactx::particles::wakefields::convolve_fft | ( | amrex::Gpu::DeviceVector< amrex::Real > const & | beam_profile_slope, |
| amrex::Gpu::DeviceVector< amrex::Real > const & | wake_func, | ||
| amrex::Real | delta_t ) |
Perform FFT-based convolution
| [in] | beam_profile_slope | number density slope along s [1/m] |
| [in] | wake_func | wake function along s [V*pc/mm] |
| [in] | delta_t | size of a bin in wake_func [m] |
| std::tuple< bool, amrex::Real > impactx::particles::wakefields::CSRBendElement | ( | VariantType const & | element_variant, |
| RefPart const & | refpart ) |
Function to calculate the radius of curvature (R) and check if an element has CSR
| [in] | element_variant | the lattice element to check |
| [in] | refpart | reference particle |
| void impactx::particles::wakefields::DepositCharge1D | ( | impactx::ImpactXParticleContainer & | myspc, |
| amrex::Gpu::DeviceVector< amrex::Real > & | charge_distribution, | ||
| amrex::Real | bin_min, | ||
| amrex::Real | bin_size, | ||
| bool | is_unity_particle_weight = false ) |
Function to calculate charge density profile
| [in] | myspc | the particle species to deposit along s |
| [out] | charge_distribution | the output array (1D) |
| [in] | bin_min | lower end of the beam in s |
| [in] | bin_size | size of the beam in s divided by num_bins |
| [in] | is_unity_particle_weight | ignore the particle weighting per macro particle, otherwise treat each particle as one physical particle |
| void impactx::particles::wakefields::DerivativeCharge1D | ( | amrex::Gpu::DeviceVector< amrex::Real > const & | charge_distribution, |
| amrex::Gpu::DeviceVector< amrex::Real > & | slopes, | ||
| amrex::Real | bin_size, | ||
| bool | GetNumberDensity = true ) |
Function to calculate the slope of the number (or charge) density
| [in] | charge_distribution | deposited charge in 1D along s |
| [out] | slopes | derivative of charge_distribution along s with len(charge_distribution)-1 |
| [in] | bin_size | size of the beam in s divided by num_bins |
| [in] | GetNumberDensity | number density if true, otherwise charge density |
| void impactx::particles::wakefields::HandleISR | ( | impactx::ImpactXParticleContainer & | particle_container, |
| T_Element const & | element_variant, | ||
| amrex::Real | slice_ds ) |
Function to handle ISR in bending dipoles.
| [in] | particle_container | the particle species container |
| [in] | element_variant | variant type of the lattice element |
| [in] | slice_ds | slice spacing along s |
| void impactx::particles::wakefields::HandleWakefield | ( | impactx::ImpactXParticleContainer & | particle_container, |
| T_Element const & | element_variant, | ||
| amrex::Real | slice_ds, | ||
| bool | print_wakefield = false ) |
Function to handle CSR bending process including charge deposition, mean transverse position calculation, wakefield generation, and convolution.
| [in] | particle_container | the particle species container |
| [in] | element_variant | variant type of the lattice element |
| [in] | slice_ds | slice spacing along s |
| [in] | print_wakefield | for debugging: print the wakefield to convolved_wakefield.txt |
| std::tuple< bool, amrex::Real > impactx::particles::wakefields::ISRBendElement | ( | VariantType const & | element_variant, |
| RefPart const & | refpart ) |
Function to calculate the radius of curvature (R) and check if an element has ISR
| [in] | element_variant | the lattice element to check |
| [in] | refpart | reference particle |
| void impactx::particles::wakefields::ISRPush | ( | ImpactXParticleContainer & | pc, |
| amrex::ParticleReal | slice_ds, | ||
| amrex::ParticleReal | rc, | ||
| int | isr_order ) |
Function to update particle momentum using a kick due to incoherent synchrotron radiation (ISR).
| [in,out] | pc | the particle species container |
| [in] | slice_ds | slice spacing [m] |
| [in] | rc | bending radius [m] |
| [in] | isr_order | order of quantum ISR effects |
| void impactx::particles::wakefields::MeanTransversePosition | ( | impactx::ImpactXParticleContainer & | myspc, |
| amrex::Gpu::DeviceVector< amrex::Real > & | mean_x, | ||
| amrex::Gpu::DeviceVector< amrex::Real > & | mean_y, | ||
| amrex::Real | bin_min, | ||
| amrex::Real | bin_size, | ||
| bool | is_unity_particle_weight = false ) |
Function to calculate the mean transverse positions (x and y) at each bin along s
| [in] | myspc | the particle species to deposit along s |
| [out] | mean_x | the output array for mean x positions (1D) |
| [out] | mean_y | the output array for mean y positions (1D) |
| [in] | bin_min | lower end of the beam in s |
| [in] | bin_size | size of the beam in s divided by num_bins |
| [in] | is_unity_particle_weight | ignore the particle weighting per macro particle, otherwise treat each particle as one physical particle |
| AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real impactx::particles::wakefields::unit_step | ( | amrex::Real | s | ) |
Step Function (Heavy-Side)
| s | value along s [m] |
| AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real impactx::particles::wakefields::w_l_csr | ( | amrex::Real | s, |
| amrex::Real | R, | ||
| amrex::Real const | bin_size ) |
CSR Wake Function
Eq. (28) in: E.L. Saldin et al., "On the coherent radiation of an electron bunch moving in an arc of a circle", NiMA Volume 398, Issues 2–3, Pages 373-394, 21 October (1997) https://doi.org/10.1016/S0168-9002(97)00822-X
| [in] | s | value along s [m] |
| [in] | R | band radius [m] |
| [in] | bin_size | longitudinal bin size [m] |
| amrex::Real impactx::particles::wakefields::w_l_rf | ( | amrex::Real | s, |
| amrex::Real | a, | ||
| amrex::Real | g, | ||
| amrex::Real | L ) |
Resistive Wall Wake Function (longitudinal)
Eq. (18) in: K. L.F. Bane, "Short-Range Dipole Wakefields in Accelerating Structures for the NLC," SLAC-PUB-9663, 2003
| [in] | s | value along s [m] |
| [in] | a | iris radius [m] |
| [in] | g | gap [m] |
| [in] | L | period length [m] |
| amrex::Real impactx::particles::wakefields::w_t_rf | ( | amrex::Real | s, |
| amrex::Real | a, | ||
| amrex::Real | g, | ||
| amrex::Real | L ) |
Resistive Wall Wake Function (transverse)
Eq. (17) in: K. L.F. Bane, "Short-Range Dipole Wakefields in Accelerating Structures for the NLC," SLAC-PUB-9663, 2003
| [in] | s | value along s [m] |
| [in] | a | iris radius [m] |
| [in] | g | gap [m] |
| [in] | L | period length [m] |
| void impactx::particles::wakefields::WakePush | ( | ImpactXParticleContainer & | pc, |
| amrex::Gpu::DeviceVector< amrex::Real > const & | convolved_wakefield, | ||
| amrex::ParticleReal | slice_ds, | ||
| amrex::Real | bin_size, | ||
| amrex::Real | bin_min ) |
Function to update particle momentum using wake force kick
| [in,out] | pc | the particle species container |
| [in] | convolved_wakefield | wake functions in [V*pc/mm] |
| [in] | slice_ds | slice spacing along s |
| [in] | bin_size | size of the beam in s divided by num_bins |
| [in] | bin_min | lower end of the beam in s |
|
constexpr |
Wake function constant [unitless]
Used in eq. (21) of: K. L.F. Bane, "Short-Range Dipole Wakefields in Accelerating Structures for the NLC," SLAC-PUB-9663, 2003
|
constexpr |
Free space impedance [Ohm].