at.physics.linear#

Coupled or non-coupled 4x4 linear motion

Functions

linopt(ring[, dp, refpts, get_chrom])

Linear analysis of a H/V coupled lattice (deprecated)

linopt2(ring, *args, **kwargs)

Linear analysis of an uncoupled lattice

linopt4(ring, *args, **kwargs)

Linear analysis of a H/V coupled lattice

linopt6(ring, *args, **kwargs)

Linear analysis of a fully coupled lattice using normal modes

avlinopt(ring[, dp, refpts])

Linear analysis of a lattice with average values

get_optics(ring[, refpts, dp, method])

Linear analysis of a fully coupled lattice

get_tune(ring, *[, method, dp, dct, df, orbit])

Computes the tunes using several available methods

get_chrom(ring, *[, method, dp, dct, df, cavpts])

Computes the chromaticities using several available methods

avlinopt(ring, dp=None, refpts=None, **kwargs)[source]#

Linear analysis of a lattice with average values

avlinopt() returns average beta, mu, dispersion over the lattice elements.

Parameters:
  • ring (Lattice) – Lattice description.

  • dp (float) – Momentum deviation.

  • refpts (Type[Element] | Element | Callable[[Element], bool] | str | None | int | Sequence[int] | bool | Sequence[bool] | RefptsCode) –

    Elements at which data is returned. It can be:

    1. an integer in the range [-len(ring), len(ring)-1] selecting the element according to python indexing rules. As a special case, len(ring) is allowed and refers to the end of the last element,

    2. an ordered list of such integers without duplicates,

    3. a numpy array of booleans of maximum length len(ring)+1, where selected elements are True.

Keyword Arguments:
  • dct (float) – Path lengthening. Defaults to None

  • df (float) – Deviation of RF frequency. Defaults to None

  • orbit (Orbit) – Avoids looking for the closed orbit if it is already known ((6,) array)

  • get_chrom (bool) – Compute chromaticities. Needs computing the tune at 2 different momentum deviations around the central one.

  • get_w (bool) – Computes chromatic amplitude functions (W) [4]. Needs to compute the optics at 2 different momentum deviations around the central one.

  • keep_lattice (bool) – Assume no lattice change since the previous tracking. Defaults to False

  • XYStep (float) – Step size. Default: DConstant.XYStep

  • DPStep (float) – Momentum step size. Default: DConstant.DPStep

  • twiss_in

    Initial conditions for transfer line optics. Record array as output by linopt6(), or dictionary. Keys:

    R or alpha, beta

    mandatory (2,) arrays

    closed_orbit

    Optional (6,) array, default 0

    dispersion

    Optional (6,) array, default 0

    If present, the attribute R will be used, otherwise the attributes alpha and beta will be used. All other attributes are ignored.

Returns:
  • elemdata – Linear optics at the points refered to by refpts, if refpts is None an empty lindata structure is returned.

  • avebeta – Average beta functions [\(\hat{\beta_x},\hat{\beta_y}\)] at refpts

  • avemu – Average phase advances [\(\hat{\mu_x},\hat{\mu_y}\)] at refpts

  • avedisp – Average dispersion [\(\hat{\eta_x}, \hat{\eta'_x}, \hat{\eta_y}, \hat{\eta'_y}\)] at refpts

  • avespos – Average s position at refpts

  • tune – [\(\nu_1,\nu_2\)], linear tunes for the two normal modes of linear motion [1]

  • chrom – [\(\xi_1,\xi_2\)], chromaticities

get_chrom(ring, *, method='linopt', dp=None, dct=None, df=None, cavpts=None, **kwargs)[source]#

Computes the chromaticities using several available methods

Parameters:
  • ring (Lattice) – Lattice description.

  • method (str) –

    'linopt' returns the tunes from the linopt6() function,

    'fft' tracks a single particle and computes the tunes with fft(),

    'interp_fft' tracks a single particle and computes the tunes with interpolated FFT.

  • dp (float) – Momentum deviation.

  • dct (float) – Path lengthening.

  • df (float) – Deviation of RF frequency.

  • cavpts (Type[Element] | Element | Callable[[Element], bool] | str | None | int | Sequence[int] | bool | Sequence[bool] | RefptsCode) – If None, look for ring.cavpts, or otherwise take all cavities.

Keyword Arguments:

DPStep (float) – Momentum step for differentiation Default: DConstant.DPStep

for the 'fft' and 'interp_fft' methods only:

Keyword Arguments:
  • nturns (int) – Number of turns. Default: 512

  • amplitude (float) – Amplitude of oscillation. Default: 1.E-6

  • remove_dc (bool) – Remove the mean of oscillation data. Default: True

  • num_harmonics (int) – Number of harmonic components to compute (before mask applied, default: 20)

  • fmin (float) – Lower tune bound. Default: 0

  • fmax (float) – Upper tune bound. Default: 1

  • hann (bool) – Turn on Hanning window. Default: False, Work only for 'fft'

Returns:

chromaticities (ndarray) – array([\(\xi_x,\xi_y\)])

get_optics(ring, refpts=None, dp=None, method=<function linopt6>, **kwargs)[source]#

Linear analysis of a fully coupled lattice

Parameters:
  • ring (Lattice) – Lattice description.

  • refpts (Type[Element] | Element | Callable[[Element], bool] | str | None | int | Sequence[int] | bool | Sequence[bool] | RefptsCode) –

    Elements at which data is returned. It can be:

    1. an integer in the range [-len(ring), len(ring)-1] selecting the element according to python indexing rules. As a special case, len(ring) is allowed and refers to the end of the last element,

    2. an ordered list of such integers without duplicates,

    3. a numpy array of booleans of maximum length len(ring)+1, where selected elements are True.

  • dp (float | None) – Momentum deviation.

  • method (Callable) –

    Method for linear optics:

    linopt2: no longitudinal motion, no H/V coupling,

    linopt4: no longitudinal motion, Sagan/Rubin 4D-analysis of coupled motion,

    linopt6 (default): with or without longitudinal motion, normal mode analysis

Keyword Arguments:
  • dct (float) – Path lengthening. Defaults to None

  • df (float) – Deviation of RF frequency. Defaults to None

  • orbit (Orbit) – Avoids looking for the closed orbit if it is already known ((6,) array)

  • get_chrom (bool) – Compute chromaticities. Needs computing the tune at 2 different momentum deviations around the central one.

  • get_w (bool) – Computes chromatic amplitude functions (W) [4]. Needs to compute the optics at 2 different momentum deviations around the central one.

  • keep_lattice (bool) – Assume no lattice change since the previous tracking. Defaults to False

  • XYStep (float) – Step size. Default: DConstant.XYStep

  • DPStep (float) – Momentum step size. Default: DConstant.DPStep

  • twiss_in

    Initial conditions for transfer line optics. Record array as output by linopt6(), or dictionary. Keys:

    R or alpha, beta

    mandatory (2,) arrays

    closed_orbit

    Optional (6,) array, default 0

    dispersion

    Optional (6,) array, default 0

    If present, the attribute R will be used, otherwise the attributes alpha and beta will be used. All other attributes are ignored.

Returns:
  • elemdata0 – Linear optics data at the entrance of the ring

  • ringdata – Lattice properties

  • elemdata – Linear optics at the points refered to by refpts, if refpts is None an empty lindata structure is returned.

Warning

The format of output record arrays depends on the selected method. See linopt2(), linopt4(), linopt6().

get_tune(ring, *, method='linopt', dp=None, dct=None, df=None, orbit=None, **kwargs)[source]#

Computes the tunes using several available methods

Parameters:
  • ring (Lattice) – Lattice description

  • method (str) –

    'linopt' returns the tunes from the linopt6() function,

    'fft' tracks a single particle and computes the tunes with fft,

    'interp_fft' tracks a single particle and computes the tunes with interpolated FFT.

  • dp (float) – Momentum deviation.

  • dct (float) – Path lengthening.

  • df (float) – Deviation of RF frequency.

  • orbit (Orbit) – Avoids looking for the closed orbit if it is already known ((6,) array)

for the 'fft' and 'interp_fft' methods only:

Keyword Arguments:
  • nturns (int) – Number of turns. Default: 512

  • amplitude (float) – Amplitude of oscillation. Default: 1.E-6

  • remove_dc (bool) – Remove the mean of oscillation data. Default: True

  • num_harmonics (int) – Number of harmonic components to compute (before mask applied, default: 20)

  • fmin (float) – Lower tune bound. Default: 0

  • fmax (float) – Upper tune bound. Default: 1

  • hann (bool) – Turn on Hanning window. Default: False. Work only for 'fft'

  • get_integer (bool) – Turn on integer tune (slower)

Returns:

tunes (ndarray) – array([\(\nu_x,\nu_y\)])

linopt2(ring, *args, **kwargs)[source]#

Linear analysis of an uncoupled lattice

Parameters:

ring (Lattice) – Lattice description.

Keyword Arguments:
  • refpts (Refpts) –

    Elements at which data is returned. It can be:

    1. an integer in the range [-len(ring), len(ring)-1] selecting the element according to python indexing rules. As a special case, len(ring) is allowed and refers to the end of the last element,

    2. an ordered list of such integers without duplicates,

    3. a numpy array of booleans of maximum length len(ring)+1, where selected elements are True.

  • dp (float) – Momentum deviation. Defaults to None

  • dct (float) – Path lengthening. Defaults to None

  • df (float) – Deviation of RF frequency. Defaults to None

  • orbit (Orbit) – Avoids looking for the closed orbit if it is already known ((6,) array)

  • get_chrom (bool) – Compute chromaticities. Needs computing the tune at 2 different momentum deviations around the central one.

  • get_w (bool) – Computes chromatic amplitude functions (W, WP) [4], and derivatives of the dispersion and twiss parameters versus dp. Needs to compute the optics at 2 different momentum deviations around the central one.

  • keep_lattice (bool) – Assume no lattice change since the previous tracking. Defaults to False

  • XYStep (float) – Step size. Default: DConstant.XYStep

  • DPStep (float) – Momentum step size. Default: DConstant.DPStep

  • twiss_in

    Initial conditions for transfer line optics. Record array as output by linopt6(), or dictionary. Keys:

    R or alpha, beta

    mandatory (2,) arrays

    closed_orbit

    Optional (6,) array, default 0

    dispersion

    Optional (6,) array, default 0

    If present, the attribute R will be used, otherwise the attributes alpha and beta will be used. All other attributes are ignored.

Returns:
  • elemdata0 – Linear optics data at the entrance of the ring

  • ringdata – Lattice properties

  • elemdata – Linear optics at the points refered to by refpts, if refpts is None an empty lindata structure is returned.

elemdata is a record array with fields:

s_pos

longitudinal position [m]

M

(4, 4) transfer matrix M from the beginning of ring to the entrance of the element [2]

closed_orbit

(6,) closed orbit vector

dispersion

(4,) dispersion vector

beta

\(\left[ \beta_x,\beta_y \right]\) vector

alpha

\(\left[ \alpha_x,\alpha_y \right]\) vector

mu

\(\left[ \mu_x,\mu_y \right]\), betatron phase (modulo \(2\pi\))

W

\(\left[ W_x,W_y \right]\) only if get_w is True: chromatic amplitude function

Wp

\(\left[ Wp_x,Wp_y \right]\) only if get_w is True: chromatic phase function

dalpha

(2,) alpha derivative vector (\(\Delta \alpha/ \delta_p\))

dbeta

(2,) beta derivative vector (\(\Delta \beta/ \delta_p\))

dmu

(2,) mu derivative vector (\(\Delta \mu/ \delta_p\))

ddispersion

(4,) dispersion derivative vector (\(\Delta D/ \delta_p\))

All values given at the entrance of each element specified in refpts. Field values can be obtained with either lindata[‘idx’] or lindata.idx

ringdata is a record array with fields:

tune

Fractional tunes

chromaticity

Chromaticities, only computed if get_chrom is True

References

[1] D.Edwards,L.Teng IEEE Trans.Nucl.Sci. NS-20, No.3, p.885-888 , 1973

[3] D.Sagan, D.Rubin Phys.Rev.Spec.Top.-Accelerators and beams, vol.2 (1999)

linopt4(ring, *args, **kwargs)[source]#

Linear analysis of a H/V coupled lattice

4D-analysis of coupled motion following Sagan/Rubin

Parameters:

ring (Lattice) – Lattice description.

Keyword Arguments:
  • refpts (Refpts) –

    Elements at which data is returned. It can be:

    1. an integer in the range [-len(ring), len(ring)-1] selecting the element according to python indexing rules. As a special case, len(ring) is allowed and refers to the end of the last element,

    2. an ordered list of such integers without duplicates,

    3. a numpy array of booleans of maximum length len(ring)+1, where selected elements are True.

  • dp (float) – Momentum deviation. Defaults to None

  • dct (float) – Path lengthening. Defaults to None

  • df (float) – Deviation of RF frequency. Defaults to None

  • orbit (Orbit) – Avoids looking for the closed orbit if it is already known ((6,) array)

  • get_chrom (bool) – Compute chromaticities. Needs computing the tune at 2 different momentum deviations around the central one.

  • get_w (bool) – Computes chromatic amplitude functions (W) [8]. Needs to compute the optics at 2 different momentum deviations around the central one.

  • keep_lattice (bool) – Assume no lattice change since the previous tracking. Defaults to False

  • XYStep (float) – Step size. Default: DConstant.XYStep

  • DPStep (float) – Momentum step size. Default: DConstant.DPStep

  • twiss_in

    Initial conditions for transfer line optics. Record array as output by linopt6(), or dictionary. Keys:

    R or alpha, beta

    mandatory (2,) arrays

    closed_orbit

    Optional (6,) array, default 0

    dispersion

    Optional (6,) array, default 0

    If present, the attribute R will be used, otherwise the attributes alpha and beta will be used. All other attributes are ignored.

Returns:
  • elemdata0 – Linear optics data at the entrance of the ring

  • ringdata – Lattice properties

  • elemdata – Linear optics at the points refered to by refpts, if refpts is None an empty lindata structure is returned.

elemdata is a record array with fields:

s_pos

longitudinal position [m]

M

(4, 4) transfer matrix M from the beginning of ring to the entrance of the element [6]

closed_orbit

(6,) closed orbit vector

dispersion

(4,) dispersion vector

beta

\(\left[ \beta_x,\beta_y \right]\) vector

alpha

\(\left[ \alpha_x,\alpha_y \right]\) vector

mu

\(\left[ \mu_x,\mu_y \right]\), betatron phase (modulo \(2\pi\))

gamma

gamma parameter of the transformation to eigenmodes [7]

W

\(\left[ W_x,W_y \right]\) only if get_w is True: chromatic amplitude function

Wp

\(\left[ Wp_x,Wp_y \right]\) only if get_w is True: chromatic phase function

dalpha

(2,) alpha derivative vector (\(\Delta \alpha/ \delta_p\))

dbeta

(2,) beta derivative vector (\(\Delta \beta/ \delta_p\))

dmu

(2,) mu derivative vector (\(\Delta \mu/ \delta_p\))

ddispersion

(4,) dispersion derivative vector (\(\Delta D/ \delta_p\))

All values given at the entrance of each element specified in refpts. Field values can be obtained with either lindata[‘idx’] or lindata.idx

ringdata is a record array with fields:

tune

Fractional tunes

chromaticity

Chromaticities, only computed if get_chrom is True

References

[5] D.Edwards,L.Teng IEEE Trans.Nucl.Sci. NS-20, No.3, p.885-888 , 1973

linopt6(ring, *args, **kwargs)[source]#

Linear analysis of a fully coupled lattice using normal modes

For circular machines, linopt6() analyses

  • the 4x4 1-turn transfer matrix if ring is 4D, or

  • the 6x6 1-turn transfer matrix if ring is 6D.

For a transfer line, The “twiss_in” intput must contain either:

  • a field R, as provided by ATLINOPT6, or

  • the fields beta and alpha, as provided by linopt and linopt6

Parameters:

ring (Lattice) – Lattice description.

Keyword Arguments:
  • refpts (Refpts) –

    Elements at which data is returned. It can be:

    1. an integer in the range [-len(ring), len(ring)-1] selecting the element according to python indexing rules. As a special case, len(ring) is allowed and refers to the end of the last element,

    2. an ordered list of such integers without duplicates,

    3. a numpy array of booleans of maximum length len(ring)+1, where selected elements are True.

  • dp (float) – Momentum deviation. Defaults to None

  • dct (float) – Path lengthening. Defaults to None

  • df (float) – Deviation of RF frequency. Defaults to None

  • orbit (Orbit) – Avoids looking for the closed orbit if it is already known ((6,) array)

  • get_chrom (bool) – Compute chromaticities. Needs computing the tune at 2 different momentum deviations around the central one.

  • get_w (bool) – Computes chromatic amplitude functions (W) [11]. Needs to compute the optics at 2 different momentum deviations around the central one.

  • keep_lattice (bool) – Assume no lattice change since the previous tracking. Defaults to False

  • XYStep (float) – Step size. Default: DConstant.XYStep

  • DPStep (float) – Momentum step size. Default: DConstant.DPStep

  • twiss_in

    Initial conditions for transfer line optics. Record array as output by linopt6(), or dictionary. Keys:

    R or alpha, beta

    mandatory (2,) arrays

    closed_orbit

    Optional (6,) array, default 0

    dispersion

    Optional (6,) array, default 0

    If present, the attribute R will be used, otherwise the attributes alpha and beta will be used. All other attributes are ignored.

  • cavpts (Refpts) – Cavity location for off-momentum tuning

Returns:
  • elemdata0 – Linear optics data at the entrance of the ring

  • ringdata – Lattice properties

  • elemdata – Linear optics at the points refered to by refpts, if refpts is None an empty lindata structure is returned.

elemdata is a record array with fields:

s_pos

longitudinal position [m]

M

(6, 6) transfer matrix M from the beginning of ring to the entrance of the element

closed_orbit

(6,) closed orbit vector

dispersion

(4,) dispersion vector

A

(6,6) A-matrix

R

(3, 6, 6) R-matrices

beta

\(\left[ \beta_x,\beta_y \right]\) vector

alpha

\(\left[ \alpha_x,\alpha_y \right]\) vector

mu

\(\left[ \mu_x,\mu_y \right]\), betatron phase (modulo \(2\pi\))

W

\(\left[ W_x,W_y \right]\) only if get_w is True: chromatic amplitude function

Wp

\(\left[ Wp_x,Wp_y \right]\) only if get_w is True: chromatic phase function

dalpha

(2,) alpha derivative vector (\(\Delta \alpha/ \delta_p\))

dbeta

(2,) beta derivative vector (\(\Delta \beta/ \delta_p\))

dmu

(2,) mu derivative vector (\(\Delta \mu/ \delta_p\))

ddispersion

(4,) dispersion derivative vector (\(\Delta D/ \delta_p\))

dR

(3, 6, 6) R derivative vector (\(\Delta R/ \delta_p\))

All values given at the entrance of each element specified in refpts. Field values can be obtained with either lindata[‘idx’] or lindata.idx

ringdata is a record array with fields:

tune

Fractional tunes

chromaticity

Chromaticities, only computed if get_chrom is True

damping_time

Damping times [s] (only if radiation is ON)

References

[9] Etienne Forest, Phys. Rev. E 58, 2481 – Published 1 August 1998

[10] Andrzej Wolski, Phys. Rev. ST Accel. Beams 9, 024001 – Published 3 February 2006