Xtrack API

Line

class xtrack.Line(elements=(), element_names=None, particle_ref=None)

Beam line object. Line.element_names contains the ordered list of beam elements, Line.element_dict is a dictionary associating to each name the corresponding beam element object.

classmethod from_dict(dct, _context=None, _buffer=None, classes=())

Build a Line from a dictionary. _context and _buffer can be optionally specified.

classmethod from_sequence(nodes=None, length=None, elements=None, sequences=None, copy_elements=False, naming_scheme='{}{}', auto_reorder=False, **kwargs)

Constructs a line from a sequence definition, inserting drift spaces as needed

Parameters
  • nodes (list of Node) – Sequence definition.

  • length – Total length (in m) of line. Determines drift behind last element.

  • elements – Dictionary with named elements, which can be refered to in the sequence definion by name.

  • sequences – Dictionary with named sub-sequences, which can be refered to in the sequence definion by name.

  • copy_elements (bool) – Whether to make copies of elements or not. By default, named elements are re-used which is memory efficient but does not allow to change parameters individually.

  • naming_scheme (str) – Naming scheme to name sub-sequences. A format string accepting two names to be joined.

  • auto_reorder (bool) – If false (default), nodes must be defined in order of increasing s coordinate, otherwise an exception is thrown. If true, nodes can be defined in any order and are re-ordered as neseccary. Useful to place additional elements inside of sub-sequences.

  • kwargs – Arguments passed to constructor of the line

Returns

An instance of Line

Examples

>>> from xtrack import Line, Node, Multipole
>>> elements = {
        'quad': Multipole(length=0.3, knl=[0, +0.50]),
        'bend': Multipole(length=0.5, knl=[np.pi / 12], hxl=[np.pi / 12]),
    }
>>> sequences = {
        'arc': [Node(1, 'quad'), Node(5, 'bend')],
    }
>>> monitor = ParticlesMonitor(...)
>>>
>>> line = Line.from_sequence([
        # direct element definition
        Node(3, xt.Multipole(...)),
        Node(7, xt.Multipole(...), name='quad1'),
        Node(1, xt.Multipole(...), name='bend1', from_='quad1'),
        ...
        # using pre-defined elements by name
        Node(13, 'quad'),
        Node(14, 'quad', name='quad3'),
        Node(2, 'bend', from_='quad3', name='bend2'),
        ....
        # using nested sequences
        Node(5, 'arc', name='section_1'),
        Node(3, monitor, from_='section_1'),
    ], length = 5, elements=elements, sequences=sequences)
classmethod from_sixinput(sixinput, classes=())
classmethod from_madx_sequence(sequence, classes=(), ignored_madtypes=(), exact_drift=False, deferred_expressions=False, install_apertures=False, apply_madx_errors=False, skip_markers=False, merge_drifts=False, merge_multipoles=False)

Build a Line from a MAD-X sequence object.

property vars
property element_refs
build_tracker(**kwargs)

Build a tracker associated for the line.

property elements
filter_elements(mask=None, exclude_types_starting_with=None)

Replace with Drifts all elements satisfying a given condition.

cycle(index_first_element=None, name_first_element=None)

Cycle the line to start from a given element.

unfreeze()
copy(_context=None, _buffer=None)
to_dict()
to_pandas()
insert_element(index=None, element=None, name=None, at_s=None, s_tol=1e-06)
append_element(element, name)
get_length()
get_s_elements(mode='upstream')

Get s position for all elements

get_s_position(at_elements=None, mode='upstream')

Get s position for given elements

remove_markers(inplace=True, keep=None)
remove_inactive_multipoles(inplace=True)
remove_zero_length_drifts(inplace=True)
merge_consecutive_drifts(inplace=True)
merge_consecutive_multipoles(inplace=True)
use_simple_quadrupoles()

Replace multipoles having only the normal quadrupolar component with quadrupole elements. The element is not replaced when synchrotron radiation is active.

use_simple_bends()

Replace multipoles having only the horizontal dipolar component with dipole elements. The element is not replaced when synchrotron radiation is active.

get_elements_of_type(types)
check_aperture()

Check that all active elements have an associated aperture.

Tracker

class xtrack.Tracker(_context=None, _buffer=None, _offset=None, line=None, sequence=None, track_kernel=None, element_classes=None, particles_class=None, skip_end_turn_actions=False, reset_s_at_end_turn=True, particles_monitor_class=None, global_xy_limit=1.0, extra_headers=(), local_particle_src=None, io_buffer=None, compile=True, use_prebuilt_kernels=True, enable_pipeline_hold=False, _element_ref_data=None)

Xsuite tracker class. It is the core of the xsuite package, allows tracking particles in a given beam line. Methods to match particle distributions and to compute twiss parameters are also available.

optimize_for_tracking(compile=True, verbose=True, keep_markers=False)

Optimize the tracker for tracking speed.

find_closed_orbit(particle_co_guess=None, particle_ref=None, co_search_settings={}, delta_zeta=0, delta0=None, zeta0=None, continue_on_closed_orbit_error=False, freeze_longitudinal=False)
compute_one_turn_matrix_finite_differences(particle_on_co, steps_r_matrix=None)
twiss(particle_ref=None, delta0=None, zeta0=None, method='6d', r_sigma=0.01, nemitt_x=1e-06, nemitt_y=1e-06, delta_disp=1e-05, delta_chrom=0.0001, particle_co_guess=None, R_matrix=None, W_matrix=None, steps_r_matrix=None, co_search_settings=None, at_elements=None, at_s=None, values_at_element_exit=False, continue_on_closed_orbit_error=False, freeze_longitudinal=False, radiation_method='full', eneloss_and_damping=False, ele_start=None, ele_stop=None, twiss_init=None, particle_on_co=None, matrix_responsiveness_tol=None, matrix_stability_tol=None, symplectify=False, reverse=False, use_full_inverse=None)
survey(X0=0, Y0=0, Z0=0, theta0=0, phi0=0, psi0=0, element0=0, reverse=False)
match(vary, targets, **kwargs)

Change a set of knobs in the beamline in order to match assigned targets. See corresponding section is the Xsuite User’s guide.

filter_elements(mask=None, exclude_types_starting_with=None)

Replace with Drifts all elements satisfying a given condition.

configure_radiation(model=None, model_beamstrahlung=None, mode='deprecated')

Configure synchrotron radiation and beamstrahlung models. Choose among: None / “mean”/ “quantum”. See corresponding section is the Xsuite User’s guide.

compensate_radiation_energy_loss(delta0=0, rtot_eneloss=1e-10, max_iter=100, **kwargs)

Compensate beam energy loss from synchrotron radiation by configuring RF cavities and Multipole elements (tapering). See corresponding section is the Xsuite User’s guide.

cycle(index_first_element=None, name_first_element=None, _buffer=None, _context=None)

Cycle the line to start from a given element.

build_particles(*args, **kwargs)

Generate a particle distribution. Equivalent to xp.Particles(tracker=tracker, …) See corresponding section is the Xsuite User’s guide.

get_backtracker(_context=None, _buffer=None, global_xy_limit='from_tracker')

Build a Tracker object that backtracks in the same line.

property particle_ref: Particles
property vars
property element_refs
property enable_pipeline_hold
get_kernel_descriptions(_context=None)
resume(session)

Resume a track session that had been placed on hold.

freeze_vars(variable_names)

Freeze assigned coordinates in tracked Particles objects.

unfreeze_vars(variable_names)

Unfreeze variables previously frozen with freeze_vars.

freeze_longitudinal(state=True)

Freeze longitudinal coordinates in tracked Particles objects. See corresponding section is the Xsuite User’s guide.

start_internal_logging_for_elements_of_type(element_type, capacity)

Start internal logging for all elements of a given type. See corresponding section is the Xsuite User’s guide.

stop_internal_logging_for_elements_of_type(element_type)

Stop internal logging for all elements of a given type. See corresponding section is the Xsuite User’s guide.

to_binary_file(path)
classmethod from_binary_file(path, particles_monitor_class=None, **kwargs) Tracker

Beam elements

Drift

class xtrack.Drift(_xobject=None, **kwargs)

Bases: BeamElement

Beam element modeling a drift section. Parameters:

  • length [m]: Length of the drift section. Default is 0.

isthick = True
behaves_like_drift = True
get_backtrack_element(_context=None, _buffer=None, _offset=None)
length

Multipole

class xtrack.Multipole(order=None, knl=None, ksl=None, **kwargs)

Bases: BeamElement

Beam element modeling a thin magnetic multipole. Parameters:

  • order [int]: Horizontal shift. Default is 0.

  • knl [m^-n, array]: Normalized integrated strength of the normal components.

  • ksl [m^-n, array]: Normalized integrated strength of the skew components.

  • hxl [rad]: Rotation angle of the reference trajectory in the horizontal plane.

  • hyl [rad]: Rotation angle of the reference trajectory in the vertical plane.

  • length [m]: Length of the originating thick multipole.

get_backtrack_element(_context=None, _buffer=None, _offset=None)
hxl
hyl
inv_factorial_order
knl
ksl
length
order
radiation_flag

Cavity

class xtrack.Cavity(_xobject=None, **kwargs)

Bases: BeamElement

Beam element modeling an RF cavity. Parameters:

  • voltage [V]: Voltage of the RF cavity. Default is 0.

  • frequency [Hz]: Frequency of the RF cavity. Default is 0.

  • lag [deg]: Phase seen by the reference particle. Default is 0.

get_backtrack_element(_context=None, _buffer=None, _offset=None)
frequency
lag
voltage

RFMultipole

class xtrack.RFMultipole(order=None, knl=None, ksl=None, pn=None, ps=None, **kwargs)

Bases: BeamElement

Beam element modeling a thin modulated multipole, with strengths dependent on the z coordinate:

kn(z) = k_n cos(2pi w tau + pn/180*pi)

ks[n](z) = k_n cos(2pi w tau + pn/180*pi)

Its parameters are:

  • order [int]: Horizontal shift. Default is 0.

  • frequency [Hz]: Frequency of the RF cavity. Default is 0.

  • knl [m^-n, array]: Normalized integrated strength of the normal components.

  • ksl [m^-n, array]: Normalized integrated strength of the skew components.

  • pn [deg, array]: Phase of the normal components.

  • ps [deg, array]: Phase of the skew components.

  • voltage [V]: Longitudinal voltage. Default is 0.

  • lag [deg]: Longitudinal phase seen by the reference particle. Default is 0.

get_backtrack_element(_context=None, _buffer=None, _offset=None)
frequency
inv_factorial_order
knl
ksl
lag
order
pn
ps
voltage

DipoleEdge

class xtrack.DipoleEdge(r21=None, r43=None, h=None, e1=None, hgap=None, fint=None, **kwargs)

Bases: BeamElement

Beam element modeling a dipole edge. Parameters:

  • h [1/m]: Curvature.

  • e1 [rad]: Face angle.

  • hgap [m]: Equivalent gap.

  • fint []: Fringe integral.

get_backtrack_element(_context=None, _buffer=None, _offset=None)
e1
fint
h
hgap
r21
r43

XYShift

class xtrack.XYShift(_xobject=None, **kwargs)

Bases: BeamElement

Beam element modeling an transverse shift of the reference system. Parameters:

  • dx [m]: Horizontal shift. Default is 0.

  • dy [m]: Vertical shift. Default is 0.

get_backtrack_element(_context=None, _buffer=None, _offset=None)
dx
dy

SRotation

class xtrack.SRotation(angle=None, cos_z=None, sin_z=None, **kwargs)

Bases: BeamElement

Beam element modeling an rotation of the reference system around the s axis. Parameters:

  • angle [deg]: Rotation angle. Default is 0.

property angle
get_backtrack_element(_context=None, _buffer=None, _offset=None)
cos_z
sin_z

LimitEllipse

class xtrack.LimitEllipse(a=None, b=None, a_squ=None, b_squ=None, **kwargs)

Bases: BeamElement

to_dict()
set_half_axes(a, b)
set_half_axes_squ(a_squ, b_squ)
get_backtrack_element(_context=None, _buffer=None, _offset=None)
a_b_squ
a_squ
b_squ

LimitRect

class xtrack.LimitRect(min_x=-10000000000.0, max_x=10000000000.0, min_y=-10000000000.0, max_y=10000000000.0, **kwargs)

Bases: BeamElement

get_backtrack_element(_context=None, _buffer=None, _offset=None)
max_x
max_y
min_x
min_y

LimitRectEllipse

class xtrack.LimitRectEllipse(max_x=10000000000.0, max_y=10000000000.0, a_squ=None, b_squ=None, a=None, b=None, **kwargs)

Bases: BeamElement

get_backtrack_element(_context=None, _buffer=None, _offset=None)
set_half_axes(a, b)
set_half_axes_squ(a_squ, b_squ)
a_b_squ
a_squ
b_squ
max_x
max_y

LimitRectPolygon

class xtrack.LimitPolygon(x_vertices=None, y_vertices=None, **kwargs)

Bases: BeamElement

x_vertices
y_vertices
x_normal
y_normal
get_backtrack_element(_context=None, _buffer=None, _offset=None)
property x_closed
property y_closed
impact_point_and_normal(x_in, y_in, z_in, x_out, y_out, z_out)
property area
property centroid
resc_fac

Monitors

class xtrack.ParticlesMonitor(_context=None, _buffer=None, _offset=None, _xobject=None, start_at_turn=None, stop_at_turn=None, n_repetitions=None, repetition_period=None, num_particles=None, particle_id_range=None, auto_to_numpy=True)

Bases: BeamElement

at_element
at_turn
beta0
charge_ratio
chi
data
delta
ebe_mode
gamma0
n_records
n_repetitions
p0c
parent_particle_id
part_id_end
part_id_start
particle_id
ptau
px
py
pzeta
repetition_period
rpp
rvv
s
start_at_turn
state
stop_at_turn
weight
x
y
zeta

Beam interaction

class xtrack.BeamInteraction(name=None, interaction_process=None, length=0, isthick=None)

Bases: object

skip_in_loss_location_refinement = True
track(particles)

BeamElement base class

class xtrack.base_element.BeamElement(_xobject=None, **kwargs)
iscollective = None
init_pipeline(pipeline_manager, name, partners_names=[])
compile_kernels(*args, **kwargs)
track(particles, increment_at_element=False)
property XoStruct
copy(_context=None, _buffer=None, _offset=None)
property extra_sources
classmethod from_dict(dct, _context=None, _buffer=None, _offset=None)
move(_context=None, _buffer=None, _offset=None)
to_dict(copy_to_cpu=True)
xoinitialize(_xobject=None, _kwargs_name_check=True, **kwargs)

Loss location refinement

class xtrack.LossLocationRefinement(tracker, backtracker=None, n_theta=None, r_max=None, dr=None, ds=None, save_refine_trackers=False, allowed_backtrack_types=())

Bases: object

refine_loss_location(particles, i_apertures=None)