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.BeamType
Beam

A 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 length nmacro that 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.
source
JuTrack.BeamMethod
Beam(r::Matrix{Float64}, energy::Float64; np::Int=size(r, 1), charge::Float64=-1.0, mass::Float64=0.51099895e6, 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.

source
JuTrack.BeamMethod

Beam(r::Matrix{Float64}; energy::Float64=1e9,np::Int=size(r, 1), charge::Float64=-1.0, mass::Float64=0.51099895e6, 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.

source
JuTrack.BeamMethod
Beam(;r::Matrix{Float64}=zeros(Float64, 1,6), energy::Float64=1e9, np::Int=size(r, 1), charge::Float64=-1.0, mass::Float64=0.51099895e6, 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.

source
JuTrack.CORRECTORType
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)
source
JuTrack.CRABCAVITYType
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)
source
JuTrack.CTPSMethod
CTPS(T::Type, TPS_Dim::Int, Max_TPS_Degree::Int)

Create a truncated power series (TPS) object with type T, dimension TPS_Dim, and maximum degree Max_TPS_Degree.

source
JuTrack.CTPSMethod
CTPS(M::CTPS{T, TPS_Dim, Max_TPS_Degree})

Copy a truncated power series (TPS) object M.

source
JuTrack.CTPSMethod
CTPS(a::T, n::Int, TPS_Dim::Int, Max_TPS_Degree::Int)

Create a truncated power series (TPS) object with type T, dimension TPS_Dim, and maximum degree Max_TPS_Degree, and set the n-th variable to a.

source
JuTrack.CTPSMethod
CTPS(a::T, TPS_Dim::Int, Max_TPS_Degree::Int)

Create a truncated power series (TPS) object with type T, dimension TPS_Dim, and maximum degree Max_TPS_Degree, and set the constant term to a.

source
JuTrack.DRIFTType
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. Not implemented yet.
  • EApertures::Array{Float64,1}: elliptical apertures. Not implemented yet.

Example:

drift = DRIFT(name="D1", len=1.0)
source
JuTrack.DRIFT_SCType
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. Not implemented yet.
  • EApertures::Array{Float64,1}: elliptical apertures. Not implemented yet.
  • 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)
source
JuTrack.EdwardsTengTwissMethod
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)=EdwardsTengTwiss(betax,betay,alphax,alphay,(1.0+alphax^2)/betax,(1.0+alphay^2)/betay,
										dx,dpx,dy,dpy,mux,muy,sin(mux),cos(mux),sin(muy),cos(muy),[R11 R12;R21 R22],mode)

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.
source
JuTrack.KOCTType
KOCT(;name::String = "OCT", len::Float64 = 0.0, k3::Float64 = 0.0, 
    PolynomA::Array{Float64,1} = zeros(Float64, 4), PolynomB::Array{Float64,1} = zeros(Float64, 4), 
    MaxOrder::Int64=3, NumIntSteps::Int64 = 10, rad::Int64=0, 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 canonical octupole element. Example:

oct = KOCT(name="O1", len=0.5, k3=0.5)
source
JuTrack.KOCT_SCType
KOCT_SC(;name::String = "OCT", len::Float64 = 0.0, k3::Float64 = 0.0, 
    PolynomA::Array{Float64,1} = zeros(Float64, 4), PolynomB::Array{Float64,1} = zeros(Float64, 4), 
    MaxOrder::Int64=3, NumIntSteps::Int64 = 10, rad::Int64=0, 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 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)
source
JuTrack.KQUADType
KQUAD(;name::String = "Quad", len::Float64 = 0.0, k1::Float64 = 0.0, 
    PolynomA::Array{Float64,1} = zeros(Float64, 4), PolynomB::Array{Float64,1} = zeros(Float64, 4), 
    MaxOrder::Int64=1, NumIntSteps::Int64 = 10, rad::Int64=0, 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 canonical quadrupole element.

Example:

quad = KQUAD(name="Q1", len=0.5, k1=0.5)
source
JuTrack.KQUAD_SCType
KQUAD_SC(;name::String = "Quad", len::Float64 = 0.0, k1::Float64 = 0.0, 
    PolynomA::Array{Float64,1} = zeros(Float64, 4), PolynomB::Array{Float64,1} = zeros(Float64, 4), 
    MaxOrder::Int64=1, NumIntSteps::Int64 = 10, rad::Int64=0, 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 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)
source
JuTrack.KSEXTType
KSEXT(;name::String = "Sext", len::Float64 = 0.0, k2::Float64 = 0.0, 
    PolynomA::Array{Float64,1} = zeros(Float64, 4), PolynomB::Array{Float64,1} = zeros(Float64, 4), 
    MaxOrder::Int64=2, NumIntSteps::Int64 = 10, rad::Int64=0, 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 canonical sextupole element. Example:

sext = KSEXT(name="S1", len=0.5, k2=0.5)
source
JuTrack.KSEXT_SCType
KSEXT_SC(;name::String = "Sext", len::Float64 = 0.0, k2::Float64 = 0.0, 
    PolynomA::Array{Float64,1} = zeros(Float64, 4), PolynomB::Array{Float64,1} = zeros(Float64, 4), 
    MaxOrder::Int64=2, NumIntSteps::Int64 = 10, rad::Int64=0, 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 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)
source
JuTrack.LongitudinalWakeType
LongitudinalWake(times::AbstractVector, wakefields::AbstractVector, wakefield::Function)

Create longitudinal wake element.

source
JuTrack.MARKERType
MARKER(;name::String = "MARKER", len::Float64 = 0.0)

A marker element. Example:

marker = MARKER(name="MARKER1")
source
JuTrack.QUADType
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)
source
JuTrack.QUAD_SCType
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)
source
JuTrack.RFCAType
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)
source
JuTrack.SBENDType
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)
source
JuTrack.SBEND_SCType
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)
source
JuTrack.SOLENOIDType
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)
source
JuTrack.StrongGaussianBeamType
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.

source
JuTrack.optics4DUCMethod
optics4DUC(bx::Float64, ax::Float64, by::Float64, ay::Float64)

Construct a 4D optics element with uncoupled optics.

Arguments

  • bx::Float64: beta function in x direction
  • ax::Float64: alpha function in x direction
  • by::Float64: beta function in y direction
  • ay::Float64: alpha function in y direction

Returns

  • optics4DUC: 4D optics element with uncoupled optics
source
JuTrack.thinMULTIPOLEType
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, 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 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])
source
JuTrack.ADcomputeRDTMethod
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.

Arguments

  • ring::Array: lattice sequence
  • index::Int: index of the element to compute the RDTs
  • changed_ids::Vector{Int}: IDs of the elements with changed parameters
  • changed_elems::Vector{Element}: elements with changed parameters
  • chromatic::Bool=true: flag to compute chromatic RDTs
  • coupling::Bool=true: flag to compute coupling RDTs
  • geometric1::Bool=true: flag to compute geometric RDTs
  • geometric2::Bool=true: flag to compute geometric RDTs
  • tuneshifts::Bool=true: flag to compute tune shifts

Returns

  • d::DrivingTerms: structure of driving terms
source
JuTrack.ADlinepass!Method
ADlinepass!(line::Vector, particles::Beam, changed_idx::Vector, changed_ele::Vector)

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: a vector of beam line elements
  • particles::Beam: a beam object
  • 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`
source
JuTrack.ADlinepass_TPSA!Method
ADlinepass_TPSA!(line::Vector, rin::Vector{CTPS{T, TPS_Dim, Max_TPS_Degree}}, changed_idx::Vector, changed_ele::Vector)

Pass 6-D 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 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`
source
JuTrack.ADringpass!Method
ADringpass!(line::Vector, particles::Beam, 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: 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`
source
JuTrack.ADtwisslineMethod
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.

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.
source
JuTrack.FMAMethod
FMA(RING, beam, nturns)

do frequency map analysis (FMA) with NAFF

Arguments

  • RING: lattice
  • beam: Beam object. Avoid zero initial coordinates.
  • nturns: number of turns

Return

  • diff_nux: difference between nux1 and nux2
  • nux1: horizontal tune 1
  • nux2: horizontal tune 2
  • nuy1: vertical tune 1
  • nuy2: vertical tune 2
source
JuTrack.Gauss3_DistMethod
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 particles
  • seed::Int=3: Random seed

Return

  • Pts1::Array{Float64,2}: 6D Gaussian distribution
source
JuTrack.RBENDMethod
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)
source
JuTrack.computeRDTMethod
computeRDT(ring, index; chromatic=false, coupling=false, geometric1=false, geometric2=false, tuneshifts=false)

Compute Hamiltonian resonance driving terms (RDTs).

Arguments

  • ring::Array: lattice sequence
  • index::Int: index of the element to compute the RDTs
  • chromatic::Bool=false: flag to compute chromatic RDTs
  • coupling::Bool=false: flag to compute coupling RDTs
  • geometric1::Bool=false: flag to compute geometric RDTs
  • geometric2::Bool=false: flag to compute geometric RDTs
  • tuneshifts::Bool=false: flag to compute tune shifts

Returns

  • d::DrivingTerms: structure of driving terms
source
JuTrack.cstMethod
cst(ctps::CTPS{T, TPS_Dim, Max_TPS_Degree})

Return the constant term of a truncated power series (TPS) object ctps.

source
JuTrack.fastfindm66Function
fastfindm66(LATTICE, dp=0.0)

Find the 6x6 transfer matrix of a lattice using numerical differentiation.

Arguments

  • LATTICE: Beam line sequence.
  • dp::Float64=0.0: Momentum deviation.

Returns

  • Matrix{Float64}: 6x6 transfer matrix.
source
JuTrack.fastfindm66_refptsMethod
fastfindm66_refpts(LATTICE, dp=0.0, refpts::Vector{Int})

Find the 6x6 transfer matrix of a lattice at specified locations using numerical differentiation.

Arguments

  • LATTICE: Beam line sequence.
  • dp::Float64=0.0: Momentum deviation.
  • refpts::Vector{Int}: Indices of elements where the transfer matrix is calculated.

Returns

  • M66_refpts: 6x6 transfer matrix at specified locations.
source
JuTrack.findelemMethod
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")
source
JuTrack.findelemMethod
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)
source
JuTrack.findindexMethod
findindex(ctps::CTPS{T, TPS_Dim, Max_TPS_Degree}, indexmap::Vector{Int})

Find the index of the indexmap in the map vector.

source
JuTrack.findm66Method
findm66(seq, dp::Float64, order::Int)

Find the 6x6 transfer matrix of a sequence using TPSA.

Arguments

  • seq: Sequence of elements.
  • dp::Float64: Momentum deviation.
  • order::Int: Order of the map.

Returns

  • Matrix{Float64}: 6x6 transfer matrix.
source
JuTrack.get_2nd_moment!Method
get_2nd_moment!(beam::Beam)

Get 2nd moment of the beam. Example:

get_2nd_moment!(beam)
println(beam.moment2nd)
source
JuTrack.get_centroid!Method
get_centroid!(beam::Beam)

Get 6-D centroid of the beam. Example:

get_centroid!(beam)
println(beam.centroid)
source
JuTrack.get_emittance!Method
get_emittance!(beam::Beam)

Get emittance of the beam. Example:

get_emittance!(beam)
println(beam.emittance)
source
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)
source
JuTrack.linepass!Method
linepass!(line::Vector, particles::Beam, refpts::Vector)

Pass particles through the line element by element. Save the particles at the reference points.

Arguments

  • line::Vector: a vector of beam line elements
  • particles::Beam: a beam object
  • refpts::Vector: a vector of reference points

Returns

  • saved_particles::Vector: a vector of saved particles
source
JuTrack.linepass!Method
linepass!(line::Vector, particles::Beam)

Pass the beam through the line element by element.

Arguments

  • line::Vector: a vector of beam line elements
  • particles::Beam: a beam object
source
JuTrack.linepass_TPSA!Method
linepass_TPSA!(line::Vector, rin::Vector{CTPS{T, TPS_Dim, Max_TPS_Degree}})

Pass 6-D TPSA coordinates through the line element by element.

Arguments

  • line::Vector: a vector of beam line elements
  • rin::Vector{CTPS{T, TPSDim, MaxTPS_Degree}}: a vector of 6-D TPSA coordinates
source
JuTrack.pass!Method
pass!(ele::DRIFT, r_in::Array{Float64,1}, num_particles::Int64, particles::Beam)

This is a function to track particles through a drift element.

Arguments

  • ele::DRIFT: a drift element
  • rin::Array{Float64,1}: 6-by-numparticles array
  • num_particles::Int64: number of particles
  • particles::Beam: beam object
source
JuTrack.periodicEdwardsTengTwissMethod
periodicEdwardsTengTwiss(seq::Vector, dp, order::Int)

Calculate the Twiss parameters for a periodic lattice.

Arguments

  • seq::Vector: 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.
source
JuTrack.plinepass!Method
plinepass!(line::Vector, particles::Beam)

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: a vector of beam line elements
  • particles::Beam: a beam object
source
JuTrack.pringpass!Method
pringpass!(line::Vector, particles::Beam, 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: a vector of beam line elements
  • particles::Beam: a beam object
  • nturn::Int: number of turns
source
JuTrack.ringpass!Method
ringpass!(line::Vector, particles::Beam, nturn::Int, save::Bool)

Pass particles through the ring for nturn turns. Save the particles at each turn.

Arguments

  • line::Vector: a vector of beam line elements
  • particles::Beam: a beam object
  • nturn::Int: number of turns
  • save::Bool: Flag
source
JuTrack.ringpass!Method
ringpass!(line::Vector, particles::Beam, nturn::Int)

Pass particles through the ring for nturn turns.

Arguments

  • line::Vector: a vector of beam line elements
  • particles::Beam: a beam object
  • nturn::Int: number of turns
source
JuTrack.ringpass_TPSA!Method
ringpass_TPSA!(line::Vector, rin::Vector{CTPS{T, TPS_Dim, Max_TPS_Degree}}, nturn::Int)

Pass 6-D 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 TPSA coordinates
  • nturn::Int: number of turns
source
JuTrack.sposMethod
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
source
JuTrack.sposMethod
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
source
JuTrack.twissPropagateMethod
twissPropagate(tin::EdwardsTengTwiss,M::Matrix{Float64})

Propagate the Twiss parameters through a matrix M.

Arguments

  • tin::EdwardsTengTwiss: Input Twiss parameters.
  • M::Matrix{Float64}: Transfer matrix.

Returns

  • EdwardsTengTwiss: Output Twiss parameters.
source
JuTrack.twisslineMethod
twissline(tin::EdwardsTengTwiss,seq::Vector, dp::Float64, order::Int, endindex::Int)

Propagate the Twiss parameters through a sequence of elements.

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.
source
JuTrack.twisslineMethod
twissline(tin::EdwardsTengTwiss,seq::Vector, 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.
source
JuTrack.twissringMethod
twissring(seq::Vector, dp::Float64, order::Int)

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.

Returns

  • twis: Twiss parameters along the ring.
source