API Reference
This section contains the detailed API reference for JuTrack.jl, including all modules, functions, and types. This reference is automatically generated from the docstrings in the package.
Modules
Below are the primary modules included in JuTrack.jl.
JuTrack.AbstractElementJuTrack.BeamJuTrack.BeamJuTrack.BeamJuTrack.BeamJuTrack.BeamJuTrack.BeamJuTrack.BeamJuTrack.CORRECTORJuTrack.CRABCAVITYJuTrack.DRIFTJuTrack.DRIFT_SCJuTrack.DRIFT_SC2P5DJuTrack.DrivingTermsJuTrack.DrivingTermsTPSADJuTrack.EdwardsTengTwissJuTrack.EdwardsTengTwissJuTrack.HOTPSAJuTrack.KOCTJuTrack.KOCT_SCJuTrack.KOCT_SC2P5DJuTrack.KQUADJuTrack.KQUAD_SCJuTrack.KQUAD_SC2P5DJuTrack.KSEXTJuTrack.KSEXT_SCJuTrack.KSEXT_SC2P5DJuTrack.LBENDJuTrack.LatticeParameterJuTrack.LongitudinalRLCWakeJuTrack.LongitudinalWakeJuTrack.MARKERJuTrack.QUADJuTrack.QUAD_SCJuTrack.QUAD_SC2P5DJuTrack.RFCAJuTrack.SBENDJuTrack.SBEND_SCJuTrack.SBEND_SC2P5DJuTrack.SOLENOIDJuTrack.StrongGaussianBeamJuTrack.WIGGLERJuTrack.optics4DUCJuTrack.thinMULTIPOLEJuTrack.ADcomputeRDTJuTrack.ADlinepass!JuTrack.ADlinepass_TPSA!JuTrack.ADringpass!JuTrack.BaxisJuTrack.Beam_GaussJuTrack.Beam_GaussJuTrack.ElementRadiationJuTrack.ElossRadiationJuTrack.Gauss3_DistJuTrack.RBENDJuTrack.WigglerRadiationJuTrack._ad_effective_lengthsJuTrack.computeDrivingTermsJuTrack.computeDrivingTermsJuTrack.computeRDTJuTrack.computeRDTJuTrack.dynamic_apertureJuTrack.elradJuTrack.fastfindm66JuTrack.fastfindm66JuTrack.fastfindm66_refptsJuTrack.fastfindm66_refptsJuTrack.find_closed_orbitJuTrack.find_closed_orbitJuTrack.find_closed_orbit_4dJuTrack.find_closed_orbit_6dJuTrack.findelemJuTrack.findelemJuTrack.findm66JuTrack.findm66JuTrack.findm66_refptsJuTrack.get_2nd_moment!JuTrack.get_centroid!JuTrack.get_emittance!JuTrack.getfocJuTrack.histogram1DinZ!JuTrack.histogram1DinZ!JuTrack.initilize_6DGaussiandist!JuTrack.linepass!JuTrack.linepass!JuTrack.linepass_TPSA!JuTrack.lu_decompositionJuTrack.lu_solveJuTrack.pass!JuTrack.periodicEdwardsTengTwissJuTrack.periodicEdwardsTengTwissJuTrack.plinepass!JuTrack.pringpass!JuTrack.ringparaJuTrack.ringpass!JuTrack.ringpass!JuTrack.ringpass_TPSA!JuTrack.sposJuTrack.sposJuTrack.symplectic_conjugate_2by2JuTrack.total_lengthJuTrack.trapzJuTrack.twissPropagateJuTrack.twissPropagateJuTrack.twiss_beamJuTrack.twisslineJuTrack.twisslineJuTrack.twisslineJuTrack.twisslineJuTrack.twissringJuTrack.twissringJuTrack.twissringJuTrack.twissringJuTrack.use_exact_driftJuTrack.wigrad
JuTrack.AbstractElement — Type
AbstractElement{T}Abstract type for all elements in the lattice.
JuTrack.Beam — Type
BeamA struct that contains the information of a beam. The struct contains the following fields:
r::Matrix{Float64}: Nx6 matrix of the 6D phase space coordinates of the beam particles.np::Int: Number of particles.nmacro::Int: Number of macro particles.energy::Float64: Energy of the beam in eV.lost_flag::Vector{Int}: Vector of lengthnmacrothat contains the lost flag of each particle.charge::Float64: Charge of the beam particles. -1.0 for electrons.mass::Float64: Mass of the beam in eV.gamma::Float64: Lorentz factor of the beam particles.beta::Float64: Velocity of the beam particles in units of the speed of light.atomnum::Float64: Atomic number of the beam particles.classrad0::Float64: Classical radiation constant. (not used)radconst::Float64: Radiation constant. (not used)T0::Float64: Revolution period of the beam particles. (not used)nturn::Int: Number of turns in the simulation. (not used)znbin::Int: Number of bins in the z direction.inzindex::Vector{Int}: Index of the z bin.zhist::Vector{Float64}: Histogram of the z direction.zhist_edges::Vector{Float64}: Edges of the z histogram.temp1::Vector{Float64}: Temporary variable.temp2::Vector{Float64}: Temporary variable.temp3::Vector{Float64}: Temporary variable.temp4::Vector{Float64}: Temporary variable.temp5::Vector{Float64}: Temporary variable.emittance::Vector{Float64}: Emittance in x, y, z directions.centroid::Vector{Float64}: Centroid in x, px, y, py, z, pz directions.moment2nd::Matrix{Float64}: 2nd momentum matrix.beamsize::Vector{Float64}: Beam size in x, px, y, py, z, pz directions.current::Float64: Beam current for space charge calculation in A.
JuTrack.Beam — Method
Beam(r::Matrix{Float64}, energy::Float64; np::Int=size(r, 1), charge::Float64=-1.0, mass::Float64=m_e, atn::Float64=1.0,
emittance::Vector{Float64}=zeros(Float64, 3), centroid::Vector{Float64}=zeros(Float64, 6), T0::Float64=0.0, znbin::Int=99,
current::Float64=0.0)Construct a Beam object with coordinates of particles and beam energy. Other parameters are optional.
JuTrack.Beam — Method
Beam(r::Matrix{Float64}; energy::Float64=1e9,np::Int=size(r, 1), charge::Float64=-1.0, mass::Float64=m_e, atn::Float64=1.0, emittance::Vector{Float64}=zeros(Float64, 3), centroid::Vector{Float64}=zeros(Float64, 6), T0::Float64=0.0, znbin::Int=99, current::Float64=0.0)
Construct a Beam object with particle coordinates. Other parameters are optional.
JuTrack.Beam — Method
Beam(;r::Matrix{Float64}=zeros(Float64, 1,6), energy::Float64=1e9, np::Int=size(r, 1), charge::Float64=-1.0, mass::Float64=m_e,
atn::Float64=1.0, emittance::Vector{Float64}=zeros(Float64, 3), centroid::Vector{Float64}=zeros(Float64, 6), T0::Float64=0.0,
znbin::Int=99, current::Float64=0.0)Construct a Beam object with default parameters. All parameters are optional.
JuTrack.Beam — Method
Beam(r_GTPS::Matrix{DTPSAD{N, T}}; energy::Float64=1e9, np::Int=size(r_GTPS, 1), charge::Float64=-1.0, mass::Float64=m_e,
atn::Float64=1.0, emittance::Vector{Float64}=zeros(Float64, 3), centroid::Vector{Float64}=zeros(Float64, 6), T0::Float64=0.0,
znbin::Int=99, current::Float64=0.0) where {N, T <: Number}Construct a Beam object with coordinates of particles in TPSA (DTPSAD type) format. Other parameters are optional.
JuTrack.Beam — Method
Beam(beam::Beam)Copy a Beam object.
JuTrack.Beam — Method
Beam(r_HOTPSA::Matrix{S}; energy::Float64=1e9, np::Int=size(r_HOTPSA, 1), charge::Float64=-1.0,
mass::Float64=m_e, atn::Float64=1.0, emittance::Vector{Float64}=zeros(Float64, 3),
centroid::Vector{Float64}=zeros(Float64, 6), T0::Float64=0.0, znbin::Int=99,
current::Float64=0.0) where {S<:AbstractHOTPSA}Construct a Beam object with coordinates of particles in the fast high-order Taylor scalar format.
JuTrack.CORRECTOR — Type
CORRECTOR(;name::String = "CORRECTOR", len::Float64 = 0.0, xkick::Float64 = 0.0, ykick::Float64 = 0.0,
T1::Array{Float64,1} = zeros(6), T2::Array{Float64,1} = zeros(6), R1::Array{Float64,2} = zeros(6,6),
R2::Array{Float64,2} = zeros(6,6))A corrector element. Example:
corrector = CORRECTOR(name="C1", len=0.5, xkick=1e-3)JuTrack.CRABCAVITY — Type
CRABCAVITY(;name::String = "CRABCAVITY", len::Float64 = 0.0, volt::Float64 = 0.0, freq::Float64 = 0.0,
phi::Float64 = 0.0, errors::Array{Float64,1} = zeros(2), energy::Float64 = 1e9)A crab cavity element. Example:
crab = CRABCAVITY(name="CRAB1", len=0.5, volt=1e6, freq=1e6)JuTrack.DRIFT — Type
DRIFT(;name::String = "DRIFT", len::Float64 = 0.0, T1::Array{Float64,1} = zeros(6),
T2::Array{Float64,1} = zeros(6), R1::Array{Float64,2} = zeros(6,6), R2::Array{Float64,2} = zeros(6,6),
RApertures::Array{Float64,1} = zeros(6), EApertures::Array{Float64,1} = zeros(6))A drift element.
Arguments
- name::String: element name
- len::Float64: element length
- T1::Array{Float64,1}: misalignment at entrance
- T2::Array{Float64,1}: misalignment at exit
- R1::Array{Float64,2}: rotation at entrance
- R2::Array{Float64,2}: rotation at exit
- RApertures::Array{Float64,1}: rectangular apertures. RApertures = [xmin, xmax, ymin, ymax, 0, 0], EApertures = [ax, ay, 0, 0, 0, 0].
- EApertures::Array{Float64,1}: elliptical apertures. RApertures = [xmin, xmax, ymin, ymax, 0, 0], EApertures = [ax, ay, 0, 0, 0, 0].
Example:
drift = DRIFT(name="D1", len=1.0)JuTrack.DRIFT_SC — Type
DRIFT_SC(;name::String = "DRIFT_SC", len::Float64 = 0.0, T1::Array{Float64,1} = zeros(6),
T2::Array{Float64,1} = zeros(6), R1::Array{Float64,2} = zeros(6,6), R2::Array{Float64,2} = zeros(6,6),
RApertures::Array{Float64,1} = zeros(6), EApertures::Array{Float64,1} = zeros(6), a::Float64 = 1.0, b::Float64 = 1.0,
Nl::Int64 = 10, Nm::Int64 = 10, Nsteps::Int64=1)A drift element with space charge.
Arguments
- name::String: element name
- len::Float64: element length
- T1::Array{Float64,1}: misalignment at entrance
- T2::Array{Float64,1}: misalignment at exit
- R1::Array{Float64,2}: rotation at entrance
- R2::Array{Float64,2}: rotation at exit
- RApertures::Array{Float64,1}: rectangular apertures. RApertures = [xmin, xmax, ymin, ymax, 0, 0], EApertures = [ax, ay, 0, 0, 0, 0].
- EApertures::Array{Float64,1}: elliptical apertures. RApertures = [xmin, xmax, ymin, ymax, 0, 0], EApertures = [ax, ay, 0, 0, 0, 0].
- a::Float64: horizontal size of the perfectly conducting pipe
- b::Float64: vertical size of the perfectly conducting pipe
- Nl::Int64: number of mode in the horizontal direction
- Nm::Int64: number of mode in the vertical direction
- Nsteps::Int64: number of steps for space charge calculation. One step represents a half-kick-half.
Example:
drift = DRIFT_SC(name="D1_SC", len=0.5, a=13e-3, b=13e-3, Nl=15, Nm=15)JuTrack.DRIFT_SC2P5D — Type
DRIFT_SC2P5D(;name::String = "DRIFT_SC2P5D", len::Float64 = 0.0,
T1::Array{Float64,1} = zeros(6), T2::Array{Float64,1} = zeros(6),
R1::Array{Float64,2} = zeros(6,6), R2::Array{Float64,2} = zeros(6,6),
RApertures::Array{Float64,1} = zeros(6), EApertures::Array{Float64,1} = zeros(6),
xsize::Int64 = 64, ysize::Int64 = 64, zsize::Int64 = 32,
pipe_radius::Float64 = 1.0, xy_ratio::Float64 = 1.0,
long_avg_n::Int64 = 3, Nsteps::Int64 = 1)A drift element with JuTrack-native 2.5-D space charge. Each step uses a half-drift, one 2.5-D RB kick over the full step length, and another half-drift.
JuTrack.DrivingTerms — Type
DrivenTermsStructure to store the driving terms for the resonance driving terms calculation.
Arguments
h21000::Vector{Float64}: Driving term h21000.h30000::Vector{Float64}: Driving term h30000.h10110::Vector{Float64}: Driving term h10110.
...
h00310::Vector{Float64}: Driving term h00310.h00400::Vector{Float64}: Driving term h00400.
JuTrack.DrivingTermsTPSAD — Type
DrivingTermsTPSAD{N}Structure to store the driving terms for the resonance driving terms calculation with DTPSAD.
Arguments
h21000::Vector{DTPSAD{N, Float64}}: Driving term h21000.h30000::Vector{DTPSAD{N, Float64}}: Driving term h30000.h10110::Vector{DTPSAD{N, Float64}}: Driving term h10110.
...
h00310::Vector{DTPSAD{N, Float64}}: Driving term h00310.h00400::Vector{DTPSAD{N, Float64}}: Driving term h00400.
JuTrack.EdwardsTengTwiss — Method
EdwardsTengTwiss(betax::Float64, betay::Float64; alphax::Float64=0.0, alphay::Float64=0.0,
dx::Float64=0.0, dy::Float64=0.0, dpx::Float64=0.0, dpy::Float64=0.0,
mux::Float64=0.0, muy::Float64=0.0,
R11::Float64=0.0, R12::Float64=0.0, R21::Float64=0.0, R22::Float64=0.0,
mode::Int=1, s::Float64=0.0)Construct a EdwardsTengTwiss object with betax and betay. All other parameters are optional.
Arguments
betax::Float64: Horizontal beta function.betay::Float64: Vertical beta function.alphax::Float64=0.0: Horizontal alpha function.alphay::Float64=0.0: Vertical alpha function.dx::Float64=0.0: Horizontal dispersion.dy::Float64=0.0: Vertical dispersion.dpx::Float64=0.0: derivative of horizontal dispersion.dpy::Float64=0.0: derivative of vertical dispersion.mux::Float64=0.0: Horizontal phase advance.muy::Float64=0.0: Vertical phase advance.R11::Float64=0.0: Matrix Element R11.R12::Float64=0.0: Matrix Element R12.R21::Float64=0.0: Matrix Element R21.R22::Float64=0.0: Matrix Element R22.mode::Int=1: mode for calculation.s::Float64=0.0: Position along the beamline.
JuTrack.EdwardsTengTwiss — Method
EdwardsTengTwiss(betax::DTPSAD{N,T}, betay::DTPSAD{N,T};
alphax::DTPSAD{N,T}=zero(DTPSAD{N,T}),
alphay::DTPSAD{N,T}=zero(DTPSAD{N,T}),
dx::DTPSAD{N,T}=zero(DTPSAD{N,T}),
dy::DTPSAD{N,T}=zero(DTPSAD{N,T}),
dpx::DTPSAD{N,T}=zero(DTPSAD{N,T}),
dpy::DTPSAD{N,T}=zero(DTPSAD{N,T}),
mux::DTPSAD{N,T}=zero(DTPSAD{N,T}),
muy::DTPSAD{N,T}=zero(DTPSAD{N,T}),
R11::DTPSAD{N,T}=zero(DTPSAD{N,T}),
R12::DTPSAD{N,T}=zero(DTPSAD{N,T}),
R21::DTPSAD{N,T}=zero(DTPSAD{N,T}),
R22::DTPSAD{N,T}=zero(DTPSAD{N,T}),
mode::Int=1, s::DTPSAD{N,T}=zero(DTPSAD{N,T})) where {N, T <: Number}Construct a EdwardsTengTwiss object with betax and betay in TPSA (DTPSAD type) format. All other parameters are optional.
Arguments
betax::DTPSAD{N,T}: Horizontal beta function.betay::DTPSAD{N,T}: Vertical beta function.alphax::DTPSAD{N,T}=zero(DTPSAD{N,T}): Horizontal alpha function.alphay::DTPSAD{N,T}=zero(DTPSAD{N,T}): Vertical alpha function.dx::DTPSAD{N,T}=zero(DTPSAD{N,T}): Horizontal dispersion.dy::DTPSAD{N,T}=zero(DTPSAD{N,T}): Vertical dispersion.dpx::DTPSAD{N,T}=zero(DTPSAD{N,T}): derivative of horizontal dispersion.dpy::DTPSAD{N,T}=zero(DTPSAD{N,T}): derivative of vertical dispersion.mux::DTPSAD{N,T}=zero(DTPSAD{N,T}): Horizontal phase advance.muy::DTPSAD{N,T}=zero(DTPSAD{N,T}): Vertical phase advance.R11::DTPSAD{N,T}=zero(DTPSAD{N,T}): Matrix Element R11.R12::DTPSAD{N,T}=zero(DTPSAD{N,T}): Matrix Element R12.R21::DTPSAD{N,T}=zero(DTPSAD{N,T}): Matrix Element R21.R22::DTPSAD{N,T}=zero(DTPSAD{N,T}): Matrix Element R22.mode::Int=1: mode for calculation.
JuTrack.HOTPSA — Type
HOTPSA{N,O,T,L} <: AbstractHOTPSAFast packed high-order Taylor scalar with N variables, truncation order O, coefficient type T, and L = binomial(N + O, O) stored coefficients.
JuTrack.KOCT — Type
KOCT(;name::String = "OCT", len::Float64 = 0.0, k0::Float64 = 0.0, k1::Float64 = 0.0, k2::Float64 = 0.0, k3::Float64 = 0.0,
PolynomA::Array{Float64,1} = zeros(Float64, 4),
MaxOrder::Int64=3, NumIntSteps::Int64 = 10, rad::Int64=0, FringeQuadEntrance::Int64 = 0,
FringeQuadExit::Int64 = 0, T1::Array{Float64,1} = zeros(Float64, 6),
T2::Array{Float64,1} = zeros(Float64, 6), R1::Array{Float64,2} = zeros(Float64, 6, 6),
R2::Array{Float64,2} = zeros(Float64, 6, 6), RApertures::Array{Float64,1} = zeros(Float64, 6),
EApertures::Array{Float64,1} = zeros(Float64, 6), KickAngle::Array{Float64,1} = zeros(Float64, 2))A canonical octupole element. Example:
oct = KOCT(name="O1", len=0.5, k3=0.5)JuTrack.KOCT_SC — Type
KOCT_SC(;name::String = "OCT", len::Float64 = 0.0, k0::Float64 = 0.0, k1::Float64 = 0.0, k2::Float64 = 0.0, k3::Float64 = 0.0,
PolynomA::Array{Float64,1} = zeros(Float64, 4),
MaxOrder::Int64=3, NumIntSteps::Int64 = 10, rad::Int64=0, FringeQuadEntrance::Int64 = 0,
FringeQuadExit::Int64 = 0, T1::Array{Float64,1} = zeros(Float64, 6),
T2::Array{Float64,1} = zeros(Float64, 6), R1::Array{Float64,2} = zeros(Float64, 6, 6),
R2::Array{Float64,2} = zeros(Float64, 6, 6), RApertures::Array{Float64,1} = zeros(Float64, 6),
EApertures::Array{Float64,1} = zeros(Float64, 6), KickAngle::Array{Float64,1} = zeros(Float64, 2),
a::Float64 = 1.0, b::Float64 = 1.0, Nl::Int64 = 10, Nm::Int64 = 10, Nsteps::Int64=1)A canonical octupole element with space charge. Example:
oct = KOCT_SC(name="O1_SC", len=0.5, k3=0.5, a=13e-3, b=13e-3, Nl=15, Nm=15)JuTrack.KOCT_SC2P5D — Type
KOCT_SC2P5D(;name::String = "OCT", len::Float64 = 0.0, k3::Float64 = 0.0,
NumIntSteps::Int64 = 10, xsize::Int64 = 64, ysize::Int64 = 64,
zsize::Int64 = 32, pipe_radius::Float64 = 1.0, xy_ratio::Float64 = 1.0,
long_avg_n::Int64 = 3, Nsteps::Int64 = 1, ...)A canonical octupole with JuTrack-native 2.5-D space charge.
JuTrack.KQUAD — Type
KQUAD(;name::String = "Quad", len::Float64 = 0.0, k1::Float64 = 0.0,
PolynomA::Array{Float64,1} = zeros(Float64, 4),
MaxOrder::Int64=1, NumIntSteps::Int64 = 10, rad::Int64=0, FringeQuadEntrance::Int64 = 0,
FringeQuadExit::Int64 = 0, T1::Array{Float64,1} = zeros(Float64, 6),
T2::Array{Float64,1} = zeros(Float64, 6), R1::Array{Float64,2} = zeros(Float64, 6, 6),
R2::Array{Float64,2} = zeros(Float64, 6, 6), RApertures::Array{Float64,1} = zeros(Float64, 6),
EApertures::Array{Float64,1} = zeros(Float64, 6), KickAngle::Array{Float64,1} = zeros(Float64, 2))A canonical quadrupole element.
Example:
quad = KQUAD(name="Q1", len=0.5, k1=0.5)JuTrack.KQUAD_SC — Type
KQUAD_SC(;name::String = "Quad", len::Float64 = 0.0, k0::Float64 = 0.0, k1::Float64 = 0.0, k2::Float64 = 0.0, k3::Float64 = 0.0,
PolynomA::Array{Float64,1} = zeros(Float64, 4),
MaxOrder::Int64=1, NumIntSteps::Int64 = 10, rad::Int64=0, FringeQuadEntrance::Int64 = 0,
FringeQuadExit::Int64 = 0, T1::Array{Float64,1} = zeros(Float64, 6),
T2::Array{Float64,1} = zeros(Float64, 6), R1::Array{Float64,2} = zeros(Float64, 6, 6),
R2::Array{Float64,2} = zeros(Float64, 6, 6), RApertures::Array{Float64,1} = zeros(Float64, 6),
EApertures::Array{Float64,1} = zeros(Float64, 6), KickAngle::Array{Float64,1} = zeros(Float64, 2),
a::Float64 = 1.0, b::Float64 = 1.0, Nl::Int64 = 10, Nm::Int64 = 10, Nsteps::Int64=1)A canonical quadrupole element with space charge. Example:
quad = KQUAD_SC(name="Q1_SC", len=0.5, k1=0.5, a=13e-3, b=13e-3, Nl=15, Nm=15)JuTrack.KQUAD_SC2P5D — Type
KQUAD_SC2P5D(;name::String = "Quad", len::Float64 = 0.0, k1::Float64 = 0.0,
NumIntSteps::Int64 = 10, xsize::Int64 = 64, ysize::Int64 = 64,
zsize::Int64 = 32, pipe_radius::Float64 = 1.0, xy_ratio::Float64 = 1.0,
long_avg_n::Int64 = 3, Nsteps::Int64 = 1, ...)A canonical quadrupole with JuTrack-native 2.5-D space charge. Each SC step tracks half of the symplectic quadrupole slice, applies one 2.5-D RB kick over the step length, and tracks the second half of the slice.
JuTrack.KSEXT — Type
KSEXT(;name::String = "Sext", len::Float64 = 0.0, k0::Float64 = 0.0, k1::Float64 = 0.0, k2::Float64 = 0.0, k3::Float64 = 0.0,
PolynomA::Array{Float64,1} = zeros(Float64, 4),
MaxOrder::Int64=2, NumIntSteps::Int64 = 10, rad::Int64=0, FringeQuadEntrance::Int64 = 0,
FringeQuadExit::Int64 = 0, T1::Array{Float64,1} = zeros(Float64, 6),
T2::Array{Float64,1} = zeros(Float64, 6), R1::Array{Float64,2} = zeros(Float64, 6, 6),
R2::Array{Float64,2} = zeros(Float64, 6, 6), RApertures::Array{Float64,1} = zeros(Float64, 6),
EApertures::Array{Float64,1} = zeros(Float64, 6), KickAngle::Array{Float64,1} = zeros(Float64, 2))A canonical sextupole element. Example:
sext = KSEXT(name="S1", len=0.5, k2=0.5)JuTrack.KSEXT_SC — Type
KSEXT_SC(;name::String = "Sext", len::Float64 = 0.0, k0::Float64 = 0.0, k1::Float64 = 0.0, k2::Float64 = 0.0, k3::Float64 = 0.0,
PolynomA::Array{Float64,1} = zeros(Float64, 4),
MaxOrder::Int64=2, NumIntSteps::Int64 = 10, rad::Int64=0, FringeQuadEntrance::Int64 = 0,
FringeQuadExit::Int64 = 0, T1::Array{Float64,1} = zeros(Float64, 6),
T2::Array{Float64,1} = zeros(Float64, 6), R1::Array{Float64,2} = zeros(Float64, 6, 6),
R2::Array{Float64,2} = zeros(Float64, 6, 6), RApertures::Array{Float64,1} = zeros(Float64, 6),
EApertures::Array{Float64,1} = zeros(Float64, 6), KickAngle::Array{Float64,1} = zeros(Float64, 2),
a::Float64 = 1.0, b::Float64 = 1.0, Nl::Int64 = 10, Nm::Int64 = 10, Nsteps::Int64=1)A canonical sextupole element with space charge. Example:
sext = KSEXT_SC(name="S1_SC", len=0.5, k2=0.5, a=13e-3, b=13e-3, Nl=15, Nm=15)JuTrack.KSEXT_SC2P5D — Type
KSEXT_SC2P5D(;name::String = "Sext", len::Float64 = 0.0, k2::Float64 = 0.0,
NumIntSteps::Int64 = 10, xsize::Int64 = 64, ysize::Int64 = 64,
zsize::Int64 = 32, pipe_radius::Float64 = 1.0, xy_ratio::Float64 = 1.0,
long_avg_n::Int64 = 3, Nsteps::Int64 = 1, ...)A canonical sextupole with JuTrack-native 2.5-D space charge.
JuTrack.LBEND — Type
LBEND(;name::String = "Bend", len::Float64 = 0.0, angle::Float64 = 0.0, e1::Float64 = 0.0, e2::Float64 = 0.0, K::Float64 = 0.0,
fint1::Float64 = 0.0, fint2::Float64 = 0.0, FullGap::Float64 = 0.0,
T1::Array{Float64,1} = zeros(Float64, 6), T2::Array{Float64,1} = zeros(Float64, 6), R1::Array{Float64,2} = zeros(Float64, 6, 6),
R2::Array{Float64,2} = zeros(Float64, 6, 6), RApertures::Array{Float64,1} = zeros(Float64, 6), EApertures::Array{Float64,1} = zeros(Float64, 6))A sector bending magnet with linear map. Example:
bend = LBEND(name="B1", len=0.5, angle=0.5)JuTrack.LatticeParameter — Type
LatticeParameter(index, field[, subindices...])
LatticeParameter(line, index, field[, subindices...])Describe one active lattice parameter without mutating the stored lattice. index is the element index in the lattice, field is the element field to replace, and optional subindices select an entry of an array field.
Pass the lattice as the first argument when the parameter will be used inside Enzyme differentiation; this stores a concrete element template and avoids mutating or rebuilding the full lattice.
Examples:
LatticeParameter(RING, id_q1, :k1) # RING[id_q1].k1
LatticeParameter(RING, id_b1, :PolynomB, 2)JuTrack.LongitudinalRLCWake — Type
LongitudinalRLCWake(;freq::Float64=1.0e9, Rshunt::Float64=1.0e6, Q0::Float64=1.0)A longitudinal RLC wake element.
JuTrack.LongitudinalWake — Type
LongitudinalWake(times::AbstractVector, wakefields::AbstractVector, wakefield::Function)Create longitudinal wake element.
Arguments
times::AbstractVector: time pointswakefields::AbstractVector: wakefield values at the time pointsfliphalf::Float64=-1.0: flip the wakefield for positrons if needed
Returns
LongitudinalWake: the created longitudinal wake element
JuTrack.MARKER — Type
MARKER(;name::String = "MARKER", len::Float64 = 0.0)A marker element. Example:
marker = MARKER(name="MARKER1")JuTrack.QUAD — Type
QUAD(;name::String = "Quad", len::Float64 = 0.0, k1::Float64 = 0.0, rad::Int64 = 0,
T1::Array{Float64,1} = zeros(6), T2::Array{Float64,1} = zeros(6), R1::Array{Float64,2} = zeros(6,6),
R2::Array{Float64,2} = zeros(6,6), RApertures::Array{Float64,1} = zeros(6), EApertures::Array{Float64,1} = zeros(6))A quadrupole element using matrix formalism. Example:
quad = QUAD(name="Q1", len=0.5, k1=1.0)JuTrack.QUAD_SC — Type
QUAD_SC(;name::String = "Quad", len::Float64 = 0.0, k1::Float64 = 0.0, rad::Int64 = 0,
T1::Array{Float64,1} = zeros(6), T2::Array{Float64,1} = zeros(6), R1::Array{Float64,2} = zeros(6,6),
R2::Array{Float64,2} = zeros(6,6), RApertures::Array{Float64,1} = zeros(6), EApertures::Array{Float64,1} = zeros(6),
a::Float64 = 1.0, b::Float64 = 1.0, Nl::Int64 = 10, Nm::Int64 = 10, Nsteps::Int64=1)A quadrupole element with space charge. Example:
quad = QUAD_SC(name="Q1_SC", len=0.5, k1=1.0, a=13e-3, b=13e-3, Nl=15, Nm=15)JuTrack.QUAD_SC2P5D — Type
QUAD_SC2P5D(;name::String = "Quad", len::Float64 = 0.0, k1::Float64 = 0.0,
xsize::Int64 = 64, ysize::Int64 = 64, zsize::Int64 = 32,
pipe_radius::Float64 = 1.0, xy_ratio::Float64 = 1.0,
long_avg_n::Int64 = 3, Nsteps::Int64 = 1, ...)A linear quadrupole element with JuTrack-native 2.5-D space charge.
JuTrack.RFCA — Type
RFCA(;name::String = "RFCA", len::Float64 = 0.0, volt::Float64 = 0.0, freq::Float64 = 0.0, h::Float64 = 1.0,
lag::Float64 = 0.0, philag::Float64 = 0.0, energy::Float64 = 0.0)A RF cavity element. Example:
rf = RFCA(name="RF1", len=0.5, volt=1e6, freq=1e6)JuTrack.SBEND — Type
SBEND(;name::String = "SBend", len::Float64 = 0.0, angle::Float64 = 0.0, e1::Float64 = 0.0, e2::Float64 = 0.0,
PolynomA::Array{Float64,1} = zeros(Float64, 4), PolynomB::Array{Float64,1} = zeros(Float64, 4),
MaxOrder::Int64=0, NumIntSteps::Int64 = 10, rad::Int64=0, fint1::Float64 = 0.0, fint2::Float64 = 0.0,
gap::Float64 = 0.0, FringeBendEntrance::Int64 = 1, FringeBendExit::Int64 = 1, FringeQuadEntrance::Int64 = 0,
FringeQuadExit::Int64 = 0, FringeIntM0::Array{Float64,1} = zeros(Float64, 5), FringeIntP0::Array{Float64,1} = zeros(Float64, 5),
T1::Array{Float64,1} = zeros(Float64, 6), T2::Array{Float64,1} = zeros(Float64, 6), R1::Array{Float64,2} = zeros(Float64, 6, 6),
R2::Array{Float64,2} = zeros(Float64, 6, 6), RApertures::Array{Float64,1} = zeros(Float64, 6), EApertures::Array{Float64,1} = zeros(Float64, 6),
KickAngle::Array{Float64,1} = zeros(Float64, 2))A sector bending magnet. Example:
bend = SBEND(name="B1", len=0.5, angle=0.5)JuTrack.SBEND_SC — Type
SBEND_SC(;name::String = "SBend", len::Float64 = 0.0, angle::Float64 = 0.0, e1::Float64 = 0.0, e2::Float64 = 0.0,
PolynomA::Array{Float64,1} = zeros(Float64, 4), PolynomB::Array{Float64,1} = zeros(Float64, 4),
MaxOrder::Int64=0, NumIntSteps::Int64 = 10, rad::Int64=0, fint1::Float64 = 0.0, fint2::Float64 = 0.0,
gap::Float64 = 0.0, FringeBendEntrance::Int64 = 1, FringeBendExit::Int64 = 1,
FringeQuadEntrance::Int64 = 0, FringeQuadExit::Int64 = 0, FringeIntM0::Array{Float64,1} = zeros(Float64, 5),
FringeIntP0::Array{Float64,1} = zeros(Float64, 5), T1::Array{Float64,1} = zeros(Float64, 6),
T2::Array{Float64,1} = zeros(Float64, 6), R1::Array{Float64,2} = zeros(Float64, 6, 6),
R2::Array{Float64,2} = zeros(Float64, 6, 6), RApertures::Array{Float64,1} = zeros(Float64, 6),
EApertures::Array{Float64,1} = zeros(Float64, 6), KickAngle::Array{Float64,1} = zeros(Float64, 2),
a::Float64 = 1.0, b::Float64 = 1.0, Nl::Int64 = 10, Nm::Int64 = 10, Nsteps::Int64=1)A sector bending magnet with space charge. Example:
bend = SBEND_SC(name="B1_SC", len=0.5, angle=0.5, a=13e-3, b=13e-3, Nl=15, Nm=15)JuTrack.SBEND_SC2P5D — Type
SBEND_SC2P5D(;name::String = "SBend", len::Float64 = 0.0, angle::Float64 = 0.0,
NumIntSteps::Int64 = 10, xsize::Int64 = 64, ysize::Int64 = 64,
zsize::Int64 = 32, pipe_radius::Float64 = 1.0, xy_ratio::Float64 = 1.0,
long_avg_n::Int64 = 3, Nsteps::Int64 = 1, ...)A sector bending magnet with JuTrack-native 2.5-D space charge.
JuTrack.SOLENOID — Type
SOLENOID(;name::String = "Solenoid", len::Float64 = 0.0, ks::Float64 = 0.0, T1::Array{Float64,1} = zeros(6),
T2::Array{Float64,1} = zeros(6), R1::Array{Float64,2} = zeros(6,6), R2::Array{Float64,2} = zeros(6,6))A solenoid element. Example:
solenoid = SOLENOID(name="S1", len=0.5, ks=1.0)JuTrack.StrongGaussianBeam — Type
StrongGaussianBeam(charge::Float64, mass::Float64, atomnum::Float64, np::Int, energy::Float64,
op::AbstractOptics4D, bs::Vector{Float64}, nz::Int)Construct a strong beam-beam element with Gaussian distribution.
Arguments
charge::Float64: charge of the particlemass::Float64: mass of the particleatomnum::Float64: atomic number of the particlenp::Int: number of particles in the beamenergy::Float64: total energy of the beamop::AbstractOptics4D: optics at the interaction pointbs::Vector{Float64}: beam size at the interaction pointnz::Int: number of slices in z direction
Returns
StrongGaussianBeam: strong beam-beam element with Gaussian distribution
JuTrack.WIGGLER — Type
WIGGLER(;name::String = "WIGGLER", len::Float64 = 0.0, lw::Float64 = 0.0, Bmax::Float64 = 0.0,
Nsteps::Int64 = 10, By::Array{Int} = [1;1;0;1;1;0], Bx::Array{Int} = Int[], energy::Float64 = 1e9,
R1::Array{Float64,2} = zeros(6,6), R2::Array{Float64,2} = zeros(6,6),
T1::Array{Float64,1} = zeros(6), T2::Array{Float64,1} = zeros(6))A wiggler element.
- arameters:
- len: total length of the wiggler (m)
- lw: period length of the wiggler (m)
- Bmax: Peak magnetic field (T)
- Nsteps: number of integration steps
- By: wiggler harmonics for horizontal wigglers. Default [1;1;0;1;1;0]
- Bx: wiggler harmonics for vertical wigglers. Default []
- energy: reference energy (eV)
Example:
wiggler = WIGGLER(name="W1", len=1.0, lw=0.1, Bmax=1.0)JuTrack.optics4DUC — Method
optics4DUC(bx::Float64, ax::Float64, by::Float64, ay::Float64)Construct a 4D optics element with uncoupled optics.
Arguments
bx::Float64: beta function in x directionax::Float64: alpha function in x directionby::Float64: beta function in y directionay::Float64: alpha function in y direction
Returns
optics4DUC: 4D optics element with uncoupled optics
JuTrack.thinMULTIPOLE — Type
thinMULTIPOLE(;name::String = "thinMULTIPOLE", len::Float64 = 0.0, PolynomA::Array{Float64,1} = zeros(Float64, 4),
PolynomB::Array{Float64,1} = zeros(Float64, 4), MaxOrder::Int64=1, NumIntSteps::Int64 = 1, rad::Int64=0,
FringeQuadEntrance::Int64 = 0, FringeQuadExit::Int64 = 0, T1::Array{Float64,1} = zeros(Float64, 6),
T2::Array{Float64,1} = zeros(Float64, 6), R1::Array{Float64,2} = zeros(Float64, 6, 6),
R2::Array{Float64,2} = zeros(Float64, 6, 6), RApertures::Array{Float64,1} = zeros(Float64, 6),
EApertures::Array{Float64,1} = zeros(Float64, 6), KickAngle::Array{Float64,1} = zeros(Float64, 2))A thin multipole element. PolynomA and PolynomB are the skew and normal components of the multipole.
Example:
multipole = thinMULTIPOLE(name="M1", len=0.5, PolynomA=[0.0, 0.0, 0.0, 0.0], PolynomB=[0.0, 0.0, 0.0, 0.0])JuTrack.ADcomputeRDT — Method
ADcomputeRDT(ring, index, changed_ids, changed_elems; chromatic=true, coupling=true, geometric1=true, geometric2=true, tuneshifts=true)Compute Hamiltonian resonance driving terms (RDTs). This function is used for auto-differentiation with Enzyme to avoid issues with mutable structs.
Arguments
ring::Array: lattice sequenceindex::Int: index of the element to compute the RDTschanged_ids::Vector{Int}: IDs of the elements with changed parameterschanged_elems::Vector{Element}: elements with changed parameterschromatic::Bool=true: flag to compute chromatic RDTscoupling::Bool=true: flag to compute coupling RDTsgeometric1::Bool=true: flag to compute geometric RDTsgeometric2::Bool=true: flag to compute geometric RDTstuneshifts::Bool=true: flag to compute tune shifts
Returns
d::DrivingTerms: structure of driving terms
JuTrack.ADlinepass! — Method
ADlinepass!(line::Vector{<:AbstractElement{Float64}}, particles::Beam{Float64},
changed_idx::Vector{Int}, changed_ele::Vector{<:AbstractElement{Float64}})Pass particles through the line element by element. The elements in the changed_idx will be replaced by the elements in changed_ele. This is a convinent function to implement automatic differentiation that avoid directly changing parameters in line`.
Arguments
- line::Vector{<:AbstractElement{Float64}}: a vector of beam line elements
- particles::Beam{Float64}: a beam object
- changed_idx::Vector{Int}: a vector of indices of the elements to be changed
- changedele::Vector{<:AbstractElement{Float64}}: a vector of elements to replace the elements in `changedidx`
JuTrack.ADlinepass_TPSA! — Method
ADlinepass_TPSA!(line::Vector{<:AbstractElement{Float64}}, rin::Vector{CTPS{T, TPS_Dim, Max_TPS_Degree}},
changed_idx::Vector, changed_ele::Vector; E0::Float64=3e9, m0::Float64=m_e) where {T, TPS_Dim, Max_TPS_Degree}Pass 6-D high-order TPSA coordinates through the line element by element. The elements in the changed_idx will be replaced by the elements in changed_ele. This is a convinent function to implement automatic differentiation that avoid directly changing parameters in line`.
Arguments
- line::Vector: a vector of beam line elements
- rin::Vector{CTPS{T, TPSDim, MaxTPS_Degree}}: a vector of 6-D high-order TPSA coordinates
- changed_idx::Vector: a vector of indices of the elements to be changed
- changedele::Vector: a vector of elements to replace the elements in `changedidx`
- E0::Float64=3e9: reference energy in eV
- m0::Float64=m_e: rest mass in eV/c^2
JuTrack.ADringpass! — Method
ADringpass!(line::Vector, particles::Beam{Float64}, nturn::Int, changed_idx::Vector, changed_ele::Vector)Pass particles through the ring for nturn turns. The elements in the changed_idx will be replaced by the elements in changed_ele. This is a convinent function to implement automatic differentiation that avoid directly changing parameters in line`.
Arguments
- line::Vector: a vector of a ring
- particles::Beam{Float64}: a beam object
- nturn::Int: number of turns
- changed_idx::Vector: a vector of indices of the elements to be changed
- changedele::Vector: a vector of elements to replace the elements in `changedidx`
JuTrack.Baxis — Method
Baxis(elem::AbstractElement, s::AbstractVector, kw::Float64, rhoinv::Float64)Compute the magnetic field on the axis of a generic wiggler.
Arguments
elem::AbstractElement: Wiggler elements::AbstractVector: Position along the wigglerkw::Float64: Wiggler wave number (2π/Lw)rhoinv::Float64: Inverse of magnetic rigidity scaled by Bmax
Returns
bx, bz: Horizontal and vertical field components (normalized by Bmax)
JuTrack.Beam_Gauss — Method
Beam_Gauss(nmacro::Int; energy::Float64=1e9, np::Int=nmacro,
betax::Float64=1.0, alphax::Float64=0.0, emitx::Float64=1e-8,
betay::Float64=1.0, alphay::Float64=0.0, emity::Float64=1e-8,
betaz::Float64=1.0, alphaz::Float64=0.0, emitz::Float64=1e-8,
charge::Float64=-1.0, mass::Float64=m_e, atn::Float64=1.0,
centroid::Vector{Float64}=zeros(Float64, 6), T0::Float64=0.0,
znbin::Int=99, current::Float64=0.0)Construct a Beam object with a 6D Gaussian distribution generated from Courant-Snyder parameters.
Arguments
nmacro::Int: Number of macro particles.energy::Float64: Beam energy in eV (default: 1e9).np::Int: Number of real particles (default: nmacro).betax, alphax, emitx: Horizontal Courant-Snyder parameters and geometric emittance.betay, alphay, emity: Vertical Courant-Snyder parameters and geometric emittance.betaz, alphaz, emitz: Longitudinal Courant-Snyder parameters and geometric emittance.charge, mass, atn: Particle species parameters.centroid: 6D centroid offset [x, px, y, py, z, pz].T0, znbin, current: Other beam parameters.
Example
beam = Beam_Gauss(10000, energy=3e9, betax=10.0, alphax=-1.5, emitx=20e-9,
betay=5.0, alphay=0.5, emity=2e-9)JuTrack.Beam_Gauss — Method
Beam_Gauss(::Type{DTPSAD}, nmacro::Int; energy::Float64=1e9, np::Int=nmacro,
betax::Float64=1.0, alphax::Float64=0.0, emitx::Float64=1e-8,
betay::Float64=1.0, alphay::Float64=0.0, emity::Float64=1e-8,
betaz::Float64=1.0, alphaz::Float64=0.0, emitz::Float64=1e-8,
charge::Float64=-1.0, mass::Float64=m_e, atn::Float64=1.0,
centroid::Vector{Float64}=zeros(Float64, 6), T0::Float64=0.0,
znbin::Int=99, current::Float64=0.0)Construct a Beam object with a 6D Gaussian distribution (DTPSAD type for automatic differentiation) generated from Courant-Snyder parameters.
Example
set_tps_dim(3)
beam = Beam_Gauss(DTPSAD, 10000, energy=3e9, betax=10.0, emitx=20e-9)JuTrack.ElementRadiation — Method
ElementRadiation(ring::Vector{<:AbstractElement{T}}, lindata::Vector{<:EdwardsTengTwiss{T}};
UseQuadrupoles::Bool=true)Compute the radiation integrals in dipoles and quadrupoles.
Analytical integration from: EVALUATION OF SYNCHROTRON RADIATION INTEGRALS R.H. Helm, M.J. Lee, P.L. Morton and M. Sands SLAC-PUB-1193, March 1973
Arguments
ring::Vector{<:AbstractElement}: Vector of lattice elementslindata::Vector{<:EdwardsTengTwiss}: Vector of Twiss parameters at element exits (length = length(ring))UseQuadrupoles::Bool=true: Include quadrupoles in radiation integrals calculation
Returns
I1, I2, I3, I4, I5, I6, Iv: Seven radiation integrals
Note
Unlike MATLAB AT, JuTrack's twissring returns Twiss parameters at element exits only. This function handles the indexing difference automatically.
Example
twiss = twissring(ring, 0.0, 0)
I1, I2, I3, I4, I5, I6, Iv = ElementRadiation(ring, twiss)JuTrack.ElossRadiation — Method
ElossRadiation(ring::Vector{<:AbstractElement{T}}, lindata::Vector{<:EdwardsTengTwiss{T}}) where TCompute the radiation integrals for energy loss elements.
This function handles special energy loss elements like SimpleQuantumDiffusion or EnergyLoss that model radiation effects without actual bending magnets.
Arguments
ring::Vector{<:AbstractElement}: Vector of lattice elementslindata::Vector{<:EdwardsTengTwiss}: Vector of Twiss parameters at element exits
Returns
I1, I2, I3, I4, I5: Five radiation integrals for energy loss elements
Note
Currently returns zeros as JuTrack does not yet have energy loss element types. This function is a placeholder for future implementation.
Example
twiss = twissring(ring, 0.0, 0)
I1e, I2e, I3e, I4e, I5e = ElossRadiation(ring, twiss)JuTrack.Gauss3_Dist — Method
Gauss3_Dist(distparam::Vector{Float64}, Npt::Int; seed::Int=3)Generate 6D Gaussian distribution.
Arguments
distparam::Vector{Float64}: Distribution parameters of the form [sigx, sigpx, muxpx, xscale, pxscale, xmu1, xmu2, sigy, sigpy, muypy, yscale, pyscale, xmu3, xmu4, sigz, sigpz, muzpz, zscale, pzscale, xmu5, xmu6]Npt::Int: Number of particlesseed::Int=3: Random seed
Return
Pts1::Array{Float64,2}: 6D Gaussian distribution
JuTrack.RBEND — Method
RBEND(;name::String = "RBend", len::Float64 = 0.0, angle::Float64 = 0.0, PolynomA::Array{Float64,1} = zeros(Float64, 4),
PolynomB::Array{Float64,1} = zeros(Float64, 4), MaxOrder::Int64=0, NumIntSteps::Int64 = 10, rad::Int64=0,
fint1::Float64 = 0.0, fint2::Float64 = 0.0, gap::Float64 = 0.0, FringeBendEntrance::Int64 = 1,
FringeBendExit::Int64 = 1, FringeQuadEntrance::Int64 = 0, FringeQuadExit::Int64 = 0,
FringeIntM0::Array{Float64,1} = zeros(Float64, 5), FringeIntP0::Array{Float64,1} = zeros(Float64, 5),
T1::Array{Float64,1} = zeros(Float64, 6), T2::Array{Float64,1} = zeros(Float64, 6),
R1::Array{Float64,2} = zeros(Float64, 6, 6), R2::Array{Float64,2} = zeros(Float64, 6, 6),
RApertures::Array{Float64,1} = zeros(Float64, 6), EApertures::Array{Float64,1} = zeros(Float64, 6),
KickAngle::Array{Float64,1} = zeros(Float64, 2))A rectangular bending magnet. Example:
bend = RBEND(name="B1", len=0.5, angle=0.5)JuTrack.WigglerRadiation — Method
WigglerRadiation(ring::Vector{<:AbstractElement{T}}, lindata::Vector{<:EdwardsTengTwiss{T}};
energy::Float64=3e9, nstep::Int=60) where TCompute the radiation integrals in wigglers.
WigglerRadiation computes the radiation integrals for all wigglers with the following approximations:
- The self-induced dispersion is neglected in I4 and I5, but is used as a lower limit for the I5 contribution
- I1, I2 are integrated analytically
- I3 is integrated analytically for a single harmonic, numerically otherwise
Arguments
ring::Vector{<:AbstractElement}: Vector of lattice elementslindata::Vector{<:EdwardsTengTwiss}: Vector of Twiss parameters at element exitsenergy::Float64=3e9: Beam energy in eVnstep::Int=60: Number of integration steps for numerical I3 calculation
Returns
I1, I2, I3, I4, I5: Five radiation integrals for wigglers
Note
This function assumes wiggler elements have the following fields:
Lw: Wiggler period lengthBmax: Maximum magnetic fieldBy: Horizontal field harmonics (matrix where each column is [_, amplitude, _, _, harmonic, phase])Bx: Vertical field harmonics (matrix where each column is [_, amplitude, _, _, harmonic, phase])
Example
twiss = twissring(ring, 0.0, 0)
I1, I2, I3, I4, I5 = WigglerRadiation(ring, twiss, energy=3e9)JuTrack._ad_effective_lengths — Method
ADtwissline(tin::EdwardsTengTwiss,seq::Vector, dp::Float64, order::Int, refpts::Vector{Int}, changed_idx::Vector, changed_ele::Vector)Propagate the Twiss parameters through a sequence of elements. Save the results at specified locations. This function is used for automatic differentiation with Enzyme to avoid access issues.
Arguments
tin::EdwardsTengTwiss: Input Twiss parameters.seq::Vector: Sequence of elements.dp::Float64: Momentum deviation.order::Int: Order of the map. 0 for finite difference, others for TPSA.refpts::Vector{Int}: Indices of elements where the Twiss parameters are calculated.changed_idx::Vector: Indices of elements with changed parameters.changed_ele::Vector: Indices of changed parameters.
Returns
Vector{EdwardsTengTwiss}: Output Twiss parameters at specified locations.
JuTrack.computeDrivingTerms — Method
computeDrivingTerms(s::Vector{Float64}, betax::Vector{Float64}, betay::Vector{Float64},
phix::Vector{Float64}, phiy::Vector{Float64}, etax::Vector{Float64}, Lista2L::Vector{Float64}, Listb2L::Vector{Float64},
Listb3L::Vector{Float64}, Listb4L::Vector{Float64}, tune::Vector{Float64}, flags::RDTflags, nPeriods::Int)Compute the resonance driving terms for a given lattice.
Arguments
s::Vector{Float64}: Vector of s positions.betax::Vector{Float64}: Vector of beta functions in x.betay::Vector{Float64}: Vector of beta functions in y.phix::Vector{Float64}: Vector of phase advances in x.phiy::Vector{Float64}: Vector of phase advances in y.etax::Vector{Float64}: Vector of eta functions in x.Lista2L::Vector{Float64}: Vector of a2L values.Listb2L::Vector{Float64}: Vector of b2L values.Listb3L::Vector{Float64}: Vector of b3L values.Listb4L::Vector{Float64}: Vector of b4L values.tune::Vector{Float64}: Vector of tune values.flags::RDTflags: Flags for the driving terms calculation.nPeriods::Int: Number of periods in the lattice.
Returns
DrivingTerms: Structure containing the computed driving terms.
JuTrack.computeDrivingTerms — Method
computeDrivingTerms(s::Vector{DTPSAD{N, T}}, betax::Vector{DTPSAD{N, Float64}}, betay::Vector{DTPSAD{N, Float64}},
phix::Vector{DTPSAD{N, Float64}}, phiy::Vector{DTPSAD{N, Float64}}, etax::Vector{DTPSAD{N, Float64}}, Lista2L::Vector{DTPSAD{N, Float64}}, Listb2L::Vector{DTPSAD{N, Float64}},
Listb3L::Vector{DTPSAD{N, Float64}}, Listb4L::Vector{DTPSAD{N, Float64}}, tune::Vector{DTPSAD{N, Float64}}, flags::RDTflags, nPeriods::Int)Compute the resonance driving terms for a given lattice using DTPSAD.
Arguments
- The same as the previous function but with DTPSAD types.
Returns
DrivingTermsTPSAD{N}: Structure containing the computed driving terms.
JuTrack.computeRDT — Method
computeRDT(ring::Vector{<:AbstractElement{Float64}}, index::Vector{Int};
chromatic=false, coupling=false, geometric1=false, geometric2=false, tuneshifts=false, E0=3e9, m0=m_e)Compute Hamiltonian resonance driving terms (RDTs)
Arguments
ring::Vector{<:AbstractElement{Float64}}: Lattice sequenceindex::Vector{Int}: Index of the element to compute the RDTschromatic::Bool: Whether to compute chromatic RDTs (default: false)coupling::Bool: Whether to compute coupling RDTs (default: false)geometric1::Bool: Whether to compute geometric RDTs of first order (default: false)geometric2::Bool: Whether to compute geometric RDTs of second order (default: false)tuneshifts::Bool: Whether to compute tune shifts (default: false)E0::Float64: Beam energy in eV (default: 3e9)m0: Particle rest mass (default: m_e)
Returns
Vector{DrivingTerms}: Vector containing the computed driving terms for each specified index
JuTrack.computeRDT — Method
computeRDT(ring::Vector{<:AbstractElement{DTPSAD{N, T}}}, index::Vector{Int};
chromatic=false, coupling=false, geometric1=false, geometric2=false, tuneshifts=false, E0=3e9, m0=m_e) where {N, T}Compute Hamiltonian resonance driving terms (RDTs)
Arguments
ring::Vector{<:AbstractElement{DTPSAD{N, T}}}: Lattice sequenceindex::Vector{Int}: Index of the element to compute the RDTschromatic::Bool: Whether to compute chromatic RDTs (default: false)coupling::Bool: Whether to compute coupling RDTs (default: false)geometric1::Bool: Whether to compute geometric RDTs of first order (default: false)geometric2::Bool: Whether to compute geometric RDTs of second order (default: false)tuneshifts::Bool: Whether to compute tune shifts (default: false)E0::Float64: Beam energy in eV (default: 3e9)m0: Particle rest mass (default: m_e)
Returns
Vector{DrivingTermsTPSAD{NVAR()}}: Vector containing the computed driving terms for each specified index
JuTrack.dynamic_aperture — Method
dynamic_aperture(RING, nturns, amp_max, amp_step, angle_steps, E, dp)Calculate the dynamic aperture of a given lattice.
Arguments
- RING::Vector{<:AbstractElement{Float64}}: a ring lattice
- nturns::Int: number of turns to track
- amp_max::Float64: maximum amplitude [m]
- amp_step::Float64: step size of amplitude scan [m]
- angle_steps::Int: number of lines. 11, 13, 19, 21, 37 etc.
- E::Float64: beam energy [eV]
- dp::Float64: momentum spread
Returns
- DA::Array{Float64,2}: boundary points of dynamic aperture
- survived_particles::Array{Float64,2}: survived particles in x and y [m]
JuTrack.elrad — Method
elrad(elem::AbstractElement{T}, dini::EdwardsTengTwiss{T}, dend::EdwardsTengTwiss{T}) where TCompute radiation integrals for a single element.
Arguments
elem::AbstractElement: Lattice elementdini::EdwardsTengTwiss: Twiss parameters at element entrancedend::EdwardsTengTwiss: Twiss parameters at element exit
Returns
- Tuple of (di1, di2, di3, di4, di5, di6, div): Radiation integral contributions
JuTrack.fastfindm66 — Function
fastfindm66(LATTICE::Vector{<:AbstractElement{Float64}}, dp=0.0; E0::Float64=3e9, m0::Float64=m_e, orb::Vector{Float64}=zeros(6))Find the 6x6 transfer matrix of a lattice using numerical differentiation.
Arguments
LATTICE: Beam line sequence.dp::Float64=0.0: Momentum deviation.E0::Float64=3e9: Reference energy in eV.m0::Float64=m_e: Particle mass.orb::Vector{Float64}=zeros(6): Initial orbit.
Returns
M66: 6x6 transfer matrix.
JuTrack.fastfindm66 — Method
fastfindm66(LATTICE::Vector{<:AbstractElement{DTPSAD{N, T}}}, dp::Float64=0.0;
E0::Float64=3e9, m0::Float64=m_e, orb::Vector{DTPSAD{N, T}}=zeros(6)) where {N, T}Find the 6x6 transfer matrix of a lattice using numerical differentiation for DTPSAD elements.
Arguments
LATTICE: Beam line sequence.dp::Float64=0.0: Momentum deviation.E0::Float64=3e9: Reference energy in eV.m0::Float64=m_e: Particle mass.orb::Vector{DTPSAD{N, T}}=zeros(6): Initial orbit.
Returns
M66: 6x6 transfer matrix.
JuTrack.fastfindm66_refpts — Method
fastfindm66_refpts(LATTICE::Vector{<:AbstractElement{Float64}}, dp::Float64, refpts::Vector{Int};
E0::Float64=3e9, m0::Float64=m_e, orb::Vector{Float64}=zeros(6))Find the 6x6 transfer matrix at specified reference points using numerical differentiation.
Arguments
LATTICE: Beam line sequence.dp::Float64: Momentum deviation.refpts::Vector{Int}: Indices of reference points.E0::Float64=3e9: Reference energy in eV.m0::Float64=m_e: Particle mass.orb::Vector{Float64}=zeros(6): Initial orbit.
Returns
M66_refpts: 6x6 transfer matrices at each reference point.
JuTrack.fastfindm66_refpts — Method
fastfindm66_refpts(LATTICE::Vector{<:AbstractElement{DTPSAD{N, T}}}, dp::Float64, refpts::Vector{Int};
E0::Float64=3e9, m0::Float64=m_e, orb::Vector{Float64}=zeros(6)) where {N, T}Find the 6x6 transfer matrix at specified reference points using numerical differentiation for DTPSAD elements.
Arguments
LATTICE: Beam line sequence.dp::Float64: Momentum deviation.refpts::Vector{Int}: Indices of reference points.E0::Float64=3e9: Reference energy in eV.m0::Float64=m_e: Particle mass.orb::Vector{Float64}=zeros(6): Initial orbit.
Returns
M66_refpts: 6x6 transfer matrices at each reference point.
JuTrack.find_closed_orbit — Function
find_closed_orbit(line::Vector{AbstractElement{Float64}}, dp::Float64=0.0; mass::Float64=m_e, energy::Float64=1e9,
guess::Vector{Float64}=zeros(Float64, 6), max_iter::Int=20, tol::Float64=1e-8)Calculates the closed orbit for a given lattice using Newton's method.
Arguments
line::Vector{AbstractElement{Float64}}: The lattice represented as a vector of elements.dp::Float64=0.0: Relative momentum deviation.mass::Float64=m_e: Mass of the particle.energy::Float64=1e9: Energy of the particle in eV.guess::Vector{Float64}=zeros(Float64, 6): Initial guess for the closed orbit.max_iter::Int=20: Maximum number of iterations.tol::Float64=1e-8: Tolerance for convergence.
Returns
x_closed::Vector{Float64}: The closed orbit coordinates.M::Matrix{Float64}: The one-turn transfer matrix at the closed orbit
JuTrack.find_closed_orbit — Method
find_closed_orbit(line::Vector{AbstractElement{DTPSAD{N, T}}}, dp::Float64=0.0; mass::Float64=m_e, energy::Float64=1e9,
guess::Vector{DTPSAD{N, T}}=zeros(DTPSAD{N, T}, 6), max_iter::Int=20, tol::Float64=1e-8) where {N, T <: Number}Calculates the closed orbit for a given lattice using Newton's method in TPSA (DTPSAD type) format.
Arguments
line::Vector{AbstractElement{DTPSAD{N, T}}}: The lattice represented as a vector of elements.dp::Float64=0.0: Relative momentum deviation.mass::Float64=m_e: Mass of the particle.energy::Float64=1e9: Energy of the particle in eV.guess::Vector{DTPSAD{N, T}}=zeros(DTPSAD{N, T}, 6): Initial guess for the closed orbit.max_iter::Int=20: Maximum number of iterations.tol::Float64=1e-8: Tolerance for convergence.
Returns
x_closed::Vector{DTPSAD{N, T}}: The closed orbitM::Matrix{DTPSAD{N, T}}: The one-turn transfer matrix at the closed orbit
JuTrack.find_closed_orbit_4d — Method
find_closed_orbit_4d(ring::Vector; dp::Float64=0.0, x0=zeros(4), energy::Float64=3.5e9, mass::Float64=m_e,
tol::Float64=1e-8, maxiter::Int=20)Find the 4-D closed orbit of a ring using autodiff. !!! Don't use this function for automatic differentiation because it already uses AD. !!! Use this function for AD will result in second order differentiation (slow or crash).
JuTrack.find_closed_orbit_6d — Method
find_closed_orbit_6d(ring::Vector; x0=zeros(6), energy::Float64=3.5e9, mass::Float64=m_e,
tol::Float64=1e-8, maxiter::Int=20)Find the 6-D closed orbit of a ring using autodiff. To find a 6-D closed orbit, the ring must have active RF cavities. !!! Don't use this function for automatic differentiation because it already uses AD. !!! Use this function for AD will result in second order differentiation (slow or crash).
JuTrack.findelem — Method
findelem(ring::Vector, field::Symbol, value)Find the index of elements with the specified field value in the ring.
Arguments
- ring::Vector: a vector of beam line elements
- field::Symbol: the field name
- value: the value of the field
Return
- ele_index::Vector{Int}: a vector of indices of the elements with the specified field value
Example
ele_index = findelem(ring, :name, "QF")JuTrack.findelem — Method
findelem(ring::Vector, type::Type)Find the index of elements with the specified type in the ring.
Arguments
- ring::Vector: a vector of beam line elements
- type::Type: the type of the element
Return
- ele_index::Vector{Int}: a vector of indices of the elements with the specified type
Example
ele_index = findelem(ring, DRIFT)JuTrack.findm66 — Method
findm66(seq::Vector{<:AbstractElement{Float64}}, dp::Float64, order::Int; E0::Float64=3e9, m0::Float64=m_e, orb::Vector{Float64}=zeros(6))Find the 6x6 transfer matrix.
Arguments
seq::Vector{<:AbstractElement{Float64}}: Sequence of elements.dp::Float64: Relative momentum deviation.order::Int: Order of the TPSA. Iforder == 0, a fast numerical method is used.E0::Float64=3e9: Reference energy in eV.m0::Float64=m_e: Particle mass.orb::Vector{Float64}=zeros(6): Initial orbit.
Returns
Matrix{Float64}: 6x6 transfer matrix.
JuTrack.findm66 — Method
findm66(seq::Vector{<:AbstractElement{DTPSAD{N, T}}}, dp::Float64, order::Int; E0::Float64=3e9, m0::Float64=m_e, orb::Vector{Float64}=zeros(6)) where {N, T}Find the 6x6 transfer matrix for a DTPSAD lattice.
Arguments
seq::Vector{<:AbstractElement{DTPSAD{N, T}}}: Sequence of elements.dp::Float64: Relative momentum deviation.order::Int: Onlyorder == 0is supported.E0::Float64=3e9: Reference energy in eV.m0::Float64=m_e: Particle mass.orb::Vector{Float64}=zeros(6): Initial orbit.
Returns
Matrix{Float64}: 6x6 transfer matrix.
JuTrack.findm66_refpts — Method
findm66_refpts(seq::Vector{<:AbstractElement{Float64}}, dp::Float64, order::Int, refpts::Vector{Int};
E0::Float64=3e9, m0::Float64=m_e, orb::Vector{Float64}=zeros(6))Find the 6x6 transfer matrix at specified reference points using high-order TPSA.
Arguments
seq::Vector{<:AbstractElement{Float64}}: Sequence of elements.dp::Float64: Relative momentum deviation.order::Int: Order of the TPSA.refpts::Vector{Int}: Indices of reference points.E0::Float64=3e9: Reference energy in eV.m0::Float64=m_e: Particle mass.orb::Vector{Float64}=zeros(6): Initial orbit.
Returns
Mat_list::Array{Float64,3}: 6x6 transfer matrices at each reference point.
JuTrack.get_2nd_moment! — Method
get_2nd_moment!(beam::Beam)Get 2nd moment of the beam. Example:
get_2nd_moment!(beam)
println(beam.moment2nd)JuTrack.get_centroid! — Method
get_centroid!(beam::Beam)Get 6-D centroid of the beam. Example:
get_centroid!(beam)
println(beam.centroid)JuTrack.get_emittance! — Method
get_emittance!(beam::Beam)Get emittance of the beam. Example:
get_emittance!(beam)
println(beam.emittance)JuTrack.getfoc — Method
getfoc(elem::AbstractElement{T}) where TGet the focusing strength (K value) from an element.
Arguments
elem::AbstractElement: Lattice element
Returns
K::Float64: Quadrupole focusing strength (K = ∂²By/∂x²)
JuTrack.histogram1DinZ! — Method
histogram1DinZ!(beam::Beam, nbins::Int64, inzindex, zhist, zhist_edges)Histogram in z.
Example:
histogram1DinZ!(beam, beam.znbin, beam.inzindex, beam.zhist, beam.zhist_edges)JuTrack.histogram1DinZ! — Method
histogram1DinZ!(beam::Beam)Histogram in z.
Example:
histogram1DinZ!(beam)JuTrack.initilize_6DGaussiandist! — Method
initilize_6DGaussiandist!(beam::Beam, optics::AbstractOptics4D, lmap::AbstractLongitudinalMap)Initialize 6D Gaussian distribution for the beam. Example:
initilize_6DGaussiandist!(beam, optics, lmap)JuTrack.linepass! — Method
linepass!(line::Vector{<:AbstractElement{Float64}}, particles::Beam{Float64}, refpts::Vector{Int})Pass particles through the line element by element. Save the particles at the reference points.
Arguments
- line::Vector{<:AbstractElement{Float64}}: a vector of beam line elements
- particles::Beam{Float64}: a beam object
- refpts::Vector{Int}: a vector of reference points
Returns
- saved_particles::Vector: a vector of saved particles at the reference points
JuTrack.linepass! — Method
linepass!(line::Vector{<:AbstractElement{Float64}}, particles::Beam{Float64})Pass the beam through the line element by element.
Arguments
- line::Vector{<:AbstractElement{Float64}}: a vector of beam line elements
- particles::Beam{Float64}: a beam object
JuTrack.linepass_TPSA! — Method
linepass_TPSA!(line::Vector{<:AbstractElement{Float64}}, rin::Vector{CTPS{T, TPS_Dim, Max_TPS_Degree}};
E0::Float64=3e9, m0::Float64=m_e) where {T, TPS_Dim, Max_TPS_Degree}Pass 6-D high-order TPSA coordinates through the line element by element.
Arguments
- line::Vector{<:AbstractElement{Float64}}: a vector of beam line elements
- rin::Vector{CTPS{T, TPSDim, MaxTPS_Degree}}: a vector of 6-D high-order TPSA coordinates
- E0::Float64=3e9: reference energy in eV
- m0::Float64=m_e: rest mass in eV/c^2
JuTrack.lu_decomposition — Method
lu_decomposition(A)Performs an in-place LU decomposition of a square matrix A.
JuTrack.lu_solve — Method
lu_solve(LU, b)Solves the system Ax = b given the LU decomposition of A.
JuTrack.pass! — Method
pass!(ele::DRIFT, r_in::Matrix{Float64}, num_particles::Int64, particles::Beam{Float64})This is a function to track particles through a drift element.
Arguments
- ele::DRIFT: a drift element
- rin::Matrix{Float64}: numparticles-by-6 matrix
- num_particles::Int64: number of particles
- particles::Beam{Float64}: beam object
JuTrack.periodicEdwardsTengTwiss — Method
periodicEdwardsTengTwiss(seq::Vector{<:AbstractElement{Float64}}, dp::Float64, order::Int)Calculate the Twiss parameters for a periodic lattice.
Arguments
seq::Vector{<:AbstractElement{Float64}}: Sequence of elements.dp::Float64: Momentum deviation.order::Int: Order of the map. 0 for finite difference, others for TPSA.
Returns
EdwardsTengTwiss: Output Twiss parameters.
JuTrack.periodicEdwardsTengTwiss — Method
periodicEdwardsTengTwiss(seq::Vector{<:AbstractElement{DTPSAD{N, T}}}, dp::Float64, order::Int) where {N, T}Calculate the Twiss parameters for a periodic lattice with DTPSAD elements.
Arguments
seq::Vector{<:AbstractElement{DTPSAD{N, T}}}: Sequence of elements.dp::Float64: Momentum deviation.order::Int: Order of the map. Only 0 is supported for DTPSAD elements.
Returns
EdwardsTengTwiss: Output Twiss parameters.
JuTrack.plinepass! — Method
plinepass!(line::Vector{<:AbstractElement{Float64}}, particles::Beam{Float64})Pass particles through the line element by element by implementing multi-threading. The number of threads is determined by the environment variable JULIA_NUM_THREADS.
Arguments
- line::Vector{<:AbstractElement{Float64}}: a vector of beam line elements
- particles::Beam{Float64}: a beam object
JuTrack.pringpass! — Method
pringpass!(line::Vector{<:AbstractElement{Float64}}, particles::Beam{Float64}, nturn::Int)Pass particles through the ring by implementing multi-threading. The number of threads is determined by the environment variable JULIA_NUM_THREADS.
Arguments
- line::Vector{<:AbstractElement{Float64}}: a vector of beam line elements
- particles::Beam{Float64}: a beam object
- nturn::Int: number of turns
JuTrack.ringpara — Method
ringpara(ring::Vector{<:AbstractElement};
energy::Float64=3e9,
Vrf::Float64=0.0,
harm::Int=1,
freq_rf::Float64=0.0,
dp::Float64=0.0,
print_summary::Bool=true)Calculate and optionally print ring parameters including equilibrium emittance, damping times, energy spread, and RF-dependent parameters.
This function computes radiation integrals from dipoles, quadrupoles, and wigglers, then calculates equilibrium beam parameters based on synchrotron radiation theory.
Arguments
ring::Vector{<:AbstractElement}: Vector of lattice elementsenergy::Float64=3e9: Beam energy in eVVrf::Float64=0.0: RF voltage per cell [V]harm::Int=1: Harmonic numberfreq_rf::Float64=0.0: RF frequency in Hzdp::Float64=0.0: Momentum deviation (for off-momentum calculations)print_summary::Bool=true: If true, print parameter summary
Returns
- NamedTuple with all calculated parameters
Example
using JuTrack
include("src/demo/SPEAR3/spear3.jl")
ring = spear3()
# Print summary
ringpara(ring, energy=3e9, Vrf=3.2e6, harm=372, freq_rf=476e6)
# Get parameters as structure
params = ringpara(ring, energy=3e9, Vrf=3e6, harm=372, freq_rf=476e6, print_summary=false)
println("Natural emittance: ", params.emittx * 1e9, " nm⋅rad")Reference
Based on AT's ringpara function. See also:
- H. Wiedemann, "Particle Accelerator Physics"
- M. Sands, "The Physics of Electron Storage Rings"
JuTrack.ringpass! — Method
ringpass!(line::Vector{<:AbstractElement{Float64}}, particles::Beam{Float64}, nturn::Int, save::Bool)Pass particles through the ring for nturn turns. Save the particles at each turn.
Arguments
- line::Vector{<:AbstractElement{Float64}}: a vector of beam line elements
- particles::Beam{Float64}: a beam object
- nturn::Int: number of turns
- save::Bool: Flag
JuTrack.ringpass! — Method
ringpass!(line::Vector{<:AbstractElement{Float64}}, particles::Beam{Float64}, nturn::Int)Pass particles through the ring for nturn turns.
Arguments
- line::Vector{<:AbstractElement{Float64}}: a vector of beam line elements
- particles::Beam{Float64}: a beam object
- nturn::Int: number of turns
JuTrack.ringpass_TPSA! — Method
ringpass_TPSA!(line::Vector{<:AbstractElement{Float64}}, rin::Vector{CTPS{T, TPS_Dim, Max_TPS_Degree}}, nturn::Int;
E0::Float64=3e9, m0::Float64=m_e) where {T, TPS_Dim, Max_TPS_Degree}Pass 6-D high-order TPSA coordinates through the ring for nturn turns.
Arguments
- line::Vector: a vector of beam line elements
- rin::Vector{CTPS{T, TPSDim, MaxTPS_Degree}}: a vector of 6-D high-order TPSA coordinates
- nturn::Int: number of turns
- E0::Float64=3e9: reference energy in eV
- m0::Float64=m_e: rest mass in eV/c^2
JuTrack.spos — Method
spos(ring::Vector, idx::Vector)Calculate the s position of the specified elements in the lattice.
Arguments
- ring::Vector: a vector of beam line elements
- idx::Vector: a vector of indices of the elements
Return
- pos::Vector{Float64}: a vector of s positions of the specified elements
JuTrack.spos — Method
spos(ring::Vector)Calculate the s position of each element in the lattice.
Arguments
- ring::Vector: a vector of beam line elements
Return
- pos::Vector{Float64}: a vector of s positions
JuTrack.symplectic_conjugate_2by2 — Method
symplectic_conjugate_2by2(M::Matrix{T}) where TCompute the symplectic conjugate of a 2x2 matrix M.
Arguments
M::Matrix{T}: The input 2x2 matrix.
Returns
Matrix{T}: The symplectic conjugate of the input matrix.
JuTrack.total_length — Method
total_length(ring::Vector)Calculate the total length of the lattice.
Arguments
- ring::Vector: a vector of beam line elements
Return
- leng::Float64: total length of the lattice
JuTrack.trapz — Method
trapz(x::AbstractVector, y::AbstractVector)Trapezoidal numerical integration.
Arguments
x::AbstractVector: Independent variable (must be sorted)y::AbstractVector: Dependent variable
Returns
- Integral of y with respect to x using trapezoidal rule
JuTrack.twissPropagate — Method
twissPropagate(tin::EdwardsTengTwiss{Float64},M::Matrix{Float64}; elem_length::Float64=0.0)Propagate the Twiss parameters through a matrix M.
Arguments
tin::EdwardsTengTwiss{Float64}: Input Twiss parameters.M::Matrix{Float64}: Transfer matrix.elem_length::Float64=0.0: Length of element to add to s position.
Returns
EdwardsTengTwiss{Float64}: Output Twiss parameters.
JuTrack.twissPropagate — Method
twissPropagate(tin::EdwardsTengTwiss{DTPSAD{N,T}},M::Matrix{DTPSAD{N,T}}; elem_length::Union{Float64,DTPSAD{N,T}}=0.0) where {N,T}Propagate the Twiss parameters through a matrix M in TPSA (DTPSAD type) format.
Arguments
tin::EdwardsTengTwiss{DTPSAD{N,T}}: Input Twiss parameters.M::Matrix{DTPSAD{N,T}}: Transfer matrix.elem_length::Union{Float64,DTPSAD{N,T}}=0.0: Length of element to add to s position.
Returns
EdwardsTengTwiss{DTPSAD{N,T}}: Output Twiss parameters.
JuTrack.twiss_beam — Method
twiss_beam(beam::Beam)Calculate twiss parameters of the beam.
JuTrack.twissline — Method
twissline(tin::EdwardsTengTwiss{Float64},seq::Vector{<:AbstractElement{Float64}}, dp::Float64, order::Int, endindex::Int)Propagate the Twiss parameters through a sequence of elements up to a specified index.
Arguments
tin::EdwardsTengTwiss: Input Twiss parameters.seq::Vector: Sequence of elements.dp::Float64: Momentum deviation.order::Int: Order of the map. 0 for finite difference, others for TPSA.endindex::Int: Index of the last element in the sequence.
Returns
EdwardsTengTwiss: Output Twiss parameters.
JuTrack.twissline — Method
twissline(tin::EdwardsTengTwiss{Float64},seq::Vector{<:AbstractElement{Float64}}, dp::Float64, order::Int, refpts::Vector{Int})Propagate the Twiss parameters through a sequence of elements. Save the results at specified locations.
Arguments
tin::EdwardsTengTwiss: Input Twiss parameters.seq::Vector: Sequence of elements.dp::Float64: Momentum deviation.order::Int: Order of the map. 0 for finite difference, others for TPSA.refpts::Vector{Int}: Indices of elements where the Twiss parameters are calculated.
Returns
Vector{EdwardsTengTwiss}: Output Twiss parameters at specified locations.
JuTrack.twissline — Method
twissline(tin::EdwardsTengTwiss{DTPSAD{N,T}},seq::Vector{<:AbstractElement{DTPSAD{N, T}}}, dp::Float64, order::Int, endindex::Int) where {N, T}Propagate the Twiss parameters through a sequence of DTPSAD elements up to a specified index
Arguments
tin::EdwardsTengTwiss: Input Twiss parameters.seq::Vector: Sequence of elements.dp::Float64: Momentum deviation.order::Int: Order of the map. Only 0 is supported for DTPSAD elements.endindex::Int: Index of the last element in the sequence.
Returns
EdwardsTengTwiss: Output Twiss parameters.
JuTrack.twissline — Method
twissline(tin::EdwardsTengTwiss{DTPSAD{N,T}},seq::Vector{<:AbstractElement{DTPSAD{N, T}}}, dp::Float64, order::Int, refpts::Vector{Int}) where {N, T}Propagate the Twiss parameters through a sequence of DTPSAD elements up to a specified index
Arguments
tin::EdwardsTengTwiss: Input Twiss parameters.seq::Vector: Sequence of elements.dp::Float64: Momentum deviation.order::Int: Order of the map. Only 0 is supported for DTPSAD elements.refpts::Vector{Int}: Indices of elements where the Twiss parameters are calculated.
Returns
EdwardsTengTwiss: Output Twiss parameters.
JuTrack.twissring — Method
twissring(seq::Vector{<:AbstractElement{Float64}}, dp::Float64, order::Int, refpts::Vector{Int};
E0::Float64=3e9, m0::Float64=m_e)Calculate the periodic Twiss parameters along the ring at specified reference points.
Arguments
seq::Vector: Sequence of elements.dp::Float64: Momentum deviation.order::Int: Order of the map. 0 for finite difference, others for TPSA.refpts::Vector{Int}: Reference points along the ring.E0::Float64: Beam energy in eV. Default is 3 GeV.m0::Float64: Particle rest mass in eV/c². Default is electron mass.
Returns
twis: Twiss parameters along the ring.
JuTrack.twissring — Method
twissring(seq::Vector{<:AbstractElement{Float64}}, dp::Float64, order::Int;
E0::Float64=3e9, m0::Float64=m_e)Calculate the periodic Twiss parameters along the ring.
Arguments
seq::Vector: Sequence of elements.dp::Float64: Momentum deviation.order::Int: Order of the map. 0 for finite difference, others for TPSA.E0::Float64: Beam energy in eV. Default is 3 GeV.m0::Float64: Particle rest mass in eV/c². Default is electron mass.
Returns
twis: Twiss parameters along the ring.
JuTrack.twissring — Method
twissring(seq::Vector{<:AbstractElement{DTPSAD{N, T}}}, dp::Float64, order::Int, refpts::Vector{Int};
E0::Float64=3e9, m0::Float64=m_e) where {N, T}Calculate the periodic Twiss parameters along the ring with DTPSAD elements at specified reference points.
Arguments
seq::Vector: Sequence of elements.dp::Float64: Momentum deviation.order::Int: Order of the map. 0 for finite difference, others for TPSA.refpts::Vector{Int}: Reference points along the ring.E0::Float64: Beam energy in eV. Default is 3 GeV.m0::Float64: Particle rest mass in eV/c². Default is electron mass.
Returns
twis: Twiss parameters along the ring.
JuTrack.twissring — Method
twissring(seq::Vector{<:AbstractElement{DTPSAD{N, T}}}, dp::Float64, order::Int;
E0::Float64=3e9, m0::Float64=m_e) where {N, T}Calculate the periodic Twiss parameters along the ring with DTPSAD elements.
Arguments
seq::Vector: Sequence of elements.dp::Float64: Momentum deviation.order::Int: Order of the map. 0 for finite difference, others for TPSA.E0::Float64: Beam energy in eV. Default is 3 GeV.m0::Float64: Particle rest mass in eV/c². Default is electron mass.
Returns
twis: Twiss parameters along the ring.
JuTrack.use_exact_drift — Method
use_exact_drift(flag)Select the longitudinal drift model:
0: legacy JuTrack linearized drift1: exact Hamiltonian drift2: corrected linearized drift with off-momentum slip
JuTrack.wigrad — Method
wigrad(elem::AbstractElement, dini::EdwardsTengTwiss, Brho::Float64, nstep::Int)Compute radiation integrals for a single wiggler element.
Arguments
elem::AbstractElement: Wiggler elementdini::EdwardsTengTwiss: Twiss parameters at wiggler entranceBrho::Float64: Magnetic rigidity (T⋅m)nstep::Int: Number of integration steps
Returns
- Tuple of (di1, di2, di3, di4, di5): Radiation integral contributions