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

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
 

Function Documentation

◆ alpha()

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

Parameters
svalue along s [m]
Returns
value of alpha(s) in Eq. (21) of Bane, 2003

◆ convolve_fft()

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

Parameters
[in]beam_profile_slopenumber density slope along s [1/m]
[in]wake_funcwake function along s [V*pc/mm]
[in]delta_tsize of a bin in wake_func [m]
Returns
FFT convolution of beam_profile_slope & wake_func (N = len(beam_profile_slope) = len(wake_func)/2)

◆ CSRBendElement()

template<typename VariantType>
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

Parameters
[in]element_variantthe lattice element to check
[in]refpartreference particle
Returns
boolean flag indicating if the element has CSR, radius of curvature (or zero)

◆ DepositCharge1D()

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

Parameters
[in]myspcthe particle species to deposit along s
[out]charge_distributionthe output array (1D)
[in]bin_minlower end of the beam in s
[in]bin_sizesize of the beam in s divided by num_bins
[in]is_unity_particle_weightignore the particle weighting per macro particle, otherwise treat each particle as one physical particle

◆ DerivativeCharge1D()

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

Parameters
[in]charge_distributiondeposited charge in 1D along s
[out]slopesderivative of charge_distribution along s with len(charge_distribution)-1
[in]bin_sizesize of the beam in s divided by num_bins
[in]GetNumberDensitynumber density if true, otherwise charge density

◆ HandleISR()

template<typename T_Element>
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.

Parameters
[in]particle_containerthe particle species container
[in]element_variantvariant type of the lattice element
[in]slice_dsslice spacing along s

◆ HandleWakefield()

template<typename T_Element>
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.

Parameters
[in]particle_containerthe particle species container
[in]element_variantvariant type of the lattice element
[in]slice_dsslice spacing along s
[in]print_wakefieldfor debugging: print the wakefield to convolved_wakefield.txt

◆ ISRBendElement()

template<typename VariantType>
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

Parameters
[in]element_variantthe lattice element to check
[in]refpartreference particle
Returns
boolean flag indicating if the element has ISR, radius of curvature (or zero)

◆ ISRPush()

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).

Parameters
[in,out]pcthe particle species container
[in]slice_dsslice spacing [m]
[in]rcbending radius [m]
[in]isr_orderorder of quantum ISR effects

◆ MeanTransversePosition()

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

Parameters
[in]myspcthe particle species to deposit along s
[out]mean_xthe output array for mean x positions (1D)
[out]mean_ythe output array for mean y positions (1D)
[in]bin_minlower end of the beam in s
[in]bin_sizesize of the beam in s divided by num_bins
[in]is_unity_particle_weightignore the particle weighting per macro particle, otherwise treat each particle as one physical particle

◆ unit_step()

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real impactx::particles::wakefields::unit_step ( amrex::Real s)

Step Function (Heavy-Side)

Parameters
svalue along s [m]
Returns
returns 1 if s>0, else 0

◆ w_l_csr()

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

Parameters
[in]svalue along s [m]
[in]Rband radius [m]
[in]bin_sizelongitudinal bin size [m]
Returns
wake functions in [V*pc/mm]

◆ w_l_rf()

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

Parameters
[in]svalue along s [m]
[in]airis radius [m]
[in]ggap [m]
[in]Lperiod length [m]
Returns
longitudinal wake field [V/m/pC/mm]

◆ w_t_rf()

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

Parameters
[in]svalue along s [m]
[in]airis radius [m]
[in]ggap [m]
[in]Lperiod length [m]
Returns
transverse wake field [V/m/pC/mm]

◆ WakePush()

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

Parameters
[in,out]pcthe particle species container
[in]convolved_wakefieldwake functions in [V*pc/mm]
[in]slice_dsslice spacing along s
[in]bin_sizesize of the beam in s divided by num_bins
[in]bin_minlower end of the beam in s

Variable Documentation

◆ alpha_1

amrex::Real impactx::particles::wakefields::alpha_1 = 0.4648
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

◆ Z0

amrex::Real impactx::particles::wakefields::Z0 = 377
constexpr

Free space impedance [Ohm].