Core structure

Simulated Commissioning

This module contains the main data structure of pySC package built up around the at.Lattice under study.

class pySC.core.simulated_commissioning.SimulatedCommissioning(ring: Lattice)[source]

The main structure of pySC, which holds all the information about lattice error sources and errors, injection settings and its errors. The class is initialized from at.Lattice.

SC.register* functions assign uncertainties SC.apply_errors functions apply those sigmas to specific elements in the lattice SC.update* functions transfer the information to the AT elements fields

Examples

>>> import at
>>> SC = SimulatedCommissioning(at.Lattice(at.Monitor('BPM')))
>>> SC.register_bpms(np.array(0), Offset=2e-6*np.ones(2))
>>> print(SC.ORD)
0
>>> print(SC.SIG.BPM[0].Offset)
>>> SC.apply_errors()
The properties of the class are:
RING:

the AT ring lattice to which the errors will be added

IDEALRING:

the ideal AT ring lattice

INJ:

Class for definition of injection properties. See Also pySC.core.classes.Injection

SIG:

Class to store the error sigma’s. This parameter is set via register_magnets, register_bpms, register_cavities, register_supports See Also pySC.core.classes.Sigmas

ORD:

Class to store the indices of registered elements in the lattice This parameter is set via register_magnets, register_bpms, register_cavities, register_supports See Also pySC.core.classes.Indices

plot:

(default=False) a boolean flag to trigger plots

apply_errors(nsigmas: float = 2)[source]

Applies errors to cavities, injection trajectory, BPMs, circumference, support structures and magnets if the corresponding uncertainties defined in SC.SIG are set. For example, for a magnet with ordinate ord every field defined in SC.SIG.Mag{ord} will be used to generate a random number using a Gaussian distribution with a cutoff (see option below) and sigma being the value of the uncertainty field. The number will be stored in the corresponding field of the lattice structure, thus SC.RING{ord}. An exeption are bending angle errors which are stored in the BendingAngleError field. See examples in the SC.register* functions for more details.

SC.apply_errors uses the uncertainties defined in SC.SIG to generate random errors and applies them to the corresponding attributes of elements in SC.RING.

Parameters:

nsigmas -- Number of sigmas at which the Gaussian distribution of errors is truncated

See also

SC.register_magnets, SC.register_support, SC.register_bpms, SC.register_cavities

get_cm_setpoints(ords: int | List[int] | ndarray, skewness: bool) ndarray[source]

Return current dipole Corrector Magnets (CM) setpoints

Reads the setpoints of the CMs specified in ords in the dimension skewness.

Parameters:
  • ords -- Array of CM ordinates in the lattice structure (ex: SC.ORD.CM[0])

  • skewness -- boolean specifying CM dimension ([False|True] -> [hor|ver])

Returns:

CM setpoints [rad]

register_bpms(ords: ndarray, **kwargs)[source]

registers BPMs specified by the ords (element indices in the lattice) in the SC class instance and initializes all required fields in the lattice elements. The ordinates of all registered BPMs are stored in SC.ORD.BPM.

Parameters:

ords -- BPM ordinates in the lattice structure.

The BPM fields in the lattice elements are:
Noise:

2 elements array of hor./ver. turn-by-turn BPM noise uncertainties (sigmas)

NoiseCO:

2 elements array of hor./ver. orbit BPM noise uncertainties (sigmas)

CalError:

2 elements array of hor./ver. BPM calibration errors uncertainties (sigmas)

Offset:

2 elements array of individual hor./ver. BPM offsets uncertainties (sigmas)

Roll:

BPM roll around z-axis w.r.t. the support structure

SumError:

Calibration error of the sum signal. The sum signal is used to determine the beam loss location with a cutoff as defined SC.INJ.beamLostAt.

Examples

Identify the ordinates of all elements named BPM and registers them as BPMs in SC:

ords = sc_tools.ords_from_regex(SC.RING,'BPM')
SC.register_bpms(ords)

Register the BPMs specified in ords in SC and set the uncertainty of the offset to 500um in both planes. A subsequent call of SC.apply_errors would generate a random BPM offset errors with sigma=500um:

SC.register_bpms(ords, Offset=500E-6*np.ones(2))

Register the BPMs specified in ords in SC and set the uncertainty of the offset to 500um in both planes. A subsequent call of SC.apply_errors would generate a random BPM offset errors with sigma=500um and a cutoff at 3 sigmas for this error:

SC.register_bpms(ords, Offset=[500E-6*np.ones(2), 3])

Register the BPMs specified in ords in SC and set the uncertainty of the offset to 500um in both planes and a calibration error of the sum signal of 20%:

SC.register_bpms(ords, Offset=500E-6*np.ones(2), SumError=0.2)

See also

bpm_reading, ords_from_regex, SC.verify_structure, SC.apply_errors, SC.register_support, SC.update_support

register_cavities(ords: ndarray, **kwargs)[source]

Register cavities specified in ords in SC by initializing all required fields in the corresponding cavity lattice elements and storing the ordinates in SC.ORD.RF.

Parameters:

ords -- Cavity ordinates in the lattice structure.

Keyword Arguments:
  • VoltageOffset -- Offset of cavity voltage wrt. to the setpoint

  • VoltageCalError -- Calibration error of cavity voltage wrt. to the setpoint

  • FrequencyOffset -- Offset of cavity frequency wrt. to the setpoint

  • FrequencyCalError -- Calibration error of cavity frequency wrt. to the setpoint

  • TimeLagOffset -- Offset of cavity phase wrt. to the setpoint

  • TimeLagCalError -- Calibration error of cavity phase wrt. to the setpoint

Examples

Identify the ordinates of all elements named ‘CAV’ and register them as cavities in SC:

ords = sc_tools.ords_from_regex(SC.RING, 'CAV')
SC.register_cavities(ords)

Register the cavities specified in ords in SC and sets the uncertainty of the frequency offset to 1kHz. A subsequent call of SC.apply_errors would generate a random frequency offset error with sigma=1kHz:

SC.register_cavities(ords, FrequencyOffset=1E3)

Register the cavities specified in ords in SC and sets the uncertainty of the frequency offset to 1kHz. A subsequent call of SC.apply_errors would generate a random frequency offset error with sigma=1kHz and a random timelag offset error (‘phase error’) with sigma=0.3m:

SC.register_cavities(ords, FrequencyOffset=1E3, TimeLagOffset=0.3)

See also

ords_from_regex, SC.verify_structure, SC.apply_errors

register_magnets(ords: ndarray, **kwargs)[source]

Registers magnets specified by ords in the SC structure and initializes all required fields in the lattice elements. The ordinates of all registered magnets are stored in SC.ORD.Magnet.

Parameters:

ords -- Magnet ordinates in the lattice structure.

Keyword Arguments:
  • CalErrorB -- Calibration error of the PolynomB fields wrt. the corresponding setpoints. Each element of the numpy.array is the error for the corresponding element in ‘PolynomB’

  • CalErrorA -- Calibration error of the PolynomA fields wrt. the corresponding setpoints. Each element of the numpy.array is the error for the corresponding element in ‘PolynomA’

  • MagnetOffset -- 3 element array of horizontal, vertical and longitudinal magnet offsets (wrt. the support structure).

  • MagnetRoll -- 3 element array [az,ax,ay] defineing magnet roll (around z-axis), pitch (roll around x-axis) and yaw (roll around y-axis); all wrt. the support structure.

  • BendingAngleError -- Error of the main bending field (corresponding uncertainty defined with BendingAngle).

  • CF -- Flag identifying the corresponding magnet as a combined function dipole/quadrupole. That implies that the bending angle depends on the quadrupole setpoint. A variation from the design value will therefore result in a bending angle error which is added to the PolynomB[0] field.

  • HCM -- Flag identifying the corresponding magnet as a horizontal corrector magnet. The corresponding value is the horizontal CM limit and stored in SC.RING[ords].CMlimit[0]. E.g. set limit to Inf.

  • VCM -- Flag identifying the corresponding magnet as a vertical corrector magnet. The corresponding value is the vertical CM limit and stored in SC.RING[ords].CMlimit[1]. E.g. set limit to Inf.

  • SkewQuad -- Flag identifying the corresponding magnet as a skew quadrupole corrector magnet. The corresponding value is the skew quadrupole limit and stored in SC.RING[ords].SkewLimit. E.g. set limit to Inf.

  • MasterOf -- Array of ordinates to which the corresponding magnet acts as master (split magnets). The magnets at ordinates ords are identified as a split magnets each with N childs as specified in the corresponding value which must be a [N x length(ords)] array. The field calculation in SC.update_magnets uses the setpoints and errors of the master magnet to calculate the child fields. The relative bending angle error of the master magnet e.g. is applied on the corresponding child bending angle appropriately. Split quadrupole magnets with different design gradients, however, can currently not be updated correctly.

If CMs or skew quadrupole correctors are specified, the ordinates are also stored in the corresponding fields SC.ORD.CM and SC.ORD.SkewQuad, respectively.

Examples

Identify the ordinates of all elements named QF and register them in SC:

ords = sc_tools.ords_from_regex(SC.RING, 'QF')
SC.register_magnets(ords)

Register the magnets specified in ords in SC and set the uncertainty of the quadrupole component to 1E-3 and 30um horizontal and vertical offset:

SC.register_magnets(ords,
                    CalErrorB=np.array([0, 1E-3]),
                    MagnetOffset=np.array([30E-6, 30E-6 0]))

Register the magnets specified in ords in SC and set the uncertainty of the quadrupole component to 1E-3, 30um horizontal and vertical offset and 100um longitudinal offset:

SC.register_magnets(ords,
                    CalErrorB=np.array([0, 1E-3]),
                    MagnetOffset=np.array([30E-6, 30E-6, 100E-6]))

Register the magnets specified in ords in SC and set the uncertainty of the roll, pitch and yaw angle to 100urad:

SC.register_magnets(ords, Roll=np.array([100E-6, 100E-6, 100E-6]))

Register the magnets specified in ords in SC and set the uncertainty of the roll, pitch and yaw angle to 100urad and 3.4 sigmas cutoff:

SC.register_magnets(ords, Roll=[100E-6*np.ones(3), 3.4])

Register split magnets. Identify the magnets named BENDa ([1xN] array masterOrds) and the magnets named BENDb and BENDc ([2xN] array childOrds) and register the masterOrds as the master magnets of the children in the corresponding columns of childOrds. The uncertainty of the bending angle is set to 1E-4:

masterOrds = sc_tools.ords_from_regex(SC.RING,'BENDa')
childOrds  = numpy.vstack((sc_tools.ords_from_regex(SC.RING,'BENDb'), sc_tools.ords_from_regex(SC.RING,'BENDc')))
SC.register_magnets(masterOrds, BendingAngle=1E-4, MasterOf=childOrds)

Register the magnets specified in ords in SC as combined function magnets and sets the uncertainty of the quadrupole component to 1E-3:

SC.register_magnets(ords, CF=1, CalErrorB=np.array([0, 1E-3]))

Register the magnets specified in ords in SC and set the uncertainty of the skew quadrupole component to 2E-3 and the uncertainty of the sextupole component to 1E-3:

SC.register_magnets(ords, CalErrorA=np.array([0, 2E-3, 0]), CalErrorB=np.array([0, 0, 1E-3]))

Register the magnets specified in ords in SC as horizontal and vertical CMs, set their dipole uncertainties to 5% and 1%, respectively and define no CM limits:

SC.register_magnets(ords, HCM=np.inf, VCM=np.inf, CalErrorB=5E-2, CalErrorA=1E-2)

Register the magnets specified in ords in SC as horizontal and vertical CMs, set their uncertainties to 5% and 1%, respectively and their limits to 1 mrad. Furthermore, set the uncertainty of the skew quadrupole component to 2E-3 and the uncertainty of the sextupole component to 1E-3:

SC.register_magnets(ords,
                    HCM=1E-3,
                    VCM=1E-3,
                    CalErrorB=np.array([5E-2, 0, 1E-3]),
                    CalErrorA=np.array([1E-2, 2E-3, 0]))

See also

ords_from_regex, SC.update_magnets, SC.verify_structure, SC.apply_errors, SC.register_support

register_supports(support_ords: ndarray, support_type: str, **kwargs)[source]

Initializes magnet support structures such as sections, plinths and girders in SC. The function input be given as name-value pairs, starting with the structure type and structure ordinates defining start-end endpoints. Optional arguments are set as the uncertainties of e.g. girder offsets in the sigma structure SC.SIG.Support.

Parameters:
  • support_ords -- [2xN] array of ordinates defining start and end locations of N registered support structures

  • support_type -- String specifying the support structure type. Valid are ‘Plinth’, ‘Girder’ or ‘Section’

  • **kwargs -- any of those listed below

Keyword Arguments:
  • Offset -- A [1x3] array defining horizontal, vertical and longitudinal offset uncertainties for the start points or [2x3] array defining horizontal, vertical and longitudinal offset uncertainties for the start end endpoints. If end points have dedicated uncertainties, SC.apply_errors applies random offset errors of both start end endpoints of the corresponding support structure, effectively tilting the support structure. If only start points have asigned uncertainties, SC.apply_errors applies to the support structure endpoints the same offset error as to the start points, resulting in a paraxial translation of the element. Only in this case dedicated ‘Roll’ uncertainties may be given which then tilt the structure around it’s center. The actual magnet or BPM offsets resulting from the support structure offsets is calculated in SC.update_support by interpolating on a straight line between girder start- and endpoints. Note that the coordinate system change due to bending magnets are ignored in this calculation. Thus, the accuracy of the result is limited if dipole magnets are involved. This may be particularly true in case of large sections and/or longitudinal offsets.

  • Roll -- [1x3] array [az,ax,ay] defining roll (around z-axis), pitch (roll around x-axis) and yaw (roll around y-axis) angle uncertainties.

Examples

Registers the girder start end endpoints defined in ords and assigns the horizontal, vertical and longitudinal girder offset uncertainties dX, dY and dZ, respectively, to the girder start points. When the support errors are applied the girder endpoints will get the same offset error as the start points, resulting in a paraxial translation of the girder:

SC.register_support(ords, "Girder", Offset=[dX, dY, dZ])

Registers the section start- end endpoints defined in ords and assigns the horizontal and vertical section offset uncertainties dX and dY, respectively, to the start points. When the support errors are applied the section endpoints will get the same offset as the start points:

SC.register_support(ords, "Section", Offset=np.array([dX, dY, 0]))

Registers the section start- end endpoints defined in ords and assigns the horizontal and vertical section offset uncertainties dX and dY, respectively, to the start points. When the support errors are applied the section endpoints will get the same offset as the start points. Also set a 4.2 sigmas cutoff:

SC.register_support(ords, "Section", Offset=[np.array([dX, dY, 0]), 4.2])

Registers the girder start end endpoints defined in ords, assigns the roll uncertainty dPhi and the horizontal and vertical girder offset uncertainties dX1 and dY1, respectively to the start points and dX2 and dY2 to the endpoints. When the support errors are applied, all girder start- and endpoints will get random offset errors and the resulting yaw and pitch angles are calculated accordingly:

SC.register_support(ords, "Girder",
                    Offset=np.array([dX1, dY1, 0; dX2, dY2, 0]),
                    Roll=np.array([dPhi, 0, 0]))

Registers the girder start end endpoints defined in ords and assigns the horizontal, vertical and longitudinal girder offset uncertainties dX, dY and dZ, respectively, and the roll, pitch and yaw angle uncertainties az, ax and ay. When the support errors are applied the girders will experience a paraxial translation according to the offsets plus the proper rotations around the three x-, y- and z-axes:

SC.register_support(ords, "Girder", Offset=np.array([dX dY dZ]), Roll=np.array([az ax ay]));

See also

ords_from_regex, SC.update_support, SC.support_offset_and_roll, plot_support, SC.apply_errors, SC.register_magnets, update_transformation

set_cavity_setpoints(ords: int | List[int] | ndarray, setpoints: float | List[float] | ndarray, param: str, method: str = 'abs')[source]

Set RF properties to setpoints

Set the setpoints of Voltage, Frequency or TimeLag as specified in “param” of the rf cavities specified in ords. If only a single setpoint is given for multiple cavities, the setpoint is applied to all cavities.

Parameters:
  • ords -- Array of cavity ordinates in the lattice structure (SC.ORD.RF)

  • setpoints -- Setpoints (array or single value for all cavities)

  • param -- String (‘Voltage’, ‘Frequency’ or ‘TimeLag’) specifying which cavity field should be set.

  • method -- ‘abs’ (default), Use absolute setpoint ‘rel’, Use relative setpoint to nominal value ‘add’, Add setpoints to current value

Examples

Sets the time lag of all cavities registered in SC to zero:

SC.set_cavity_setpoints(ords=SC.ORD.RF, setpoints=0.0, param='TimeLag')

Adds 1kHz to the frequency of the first cavity:

SC.set_cavity_setpoints(ords=SC.ORD.RF[0], setpoints=1E3, param='Frequency', method='add')
set_cm_setpoints(ords: int | List[int] | ndarray, setpoints: float | List[float] | ndarray, skewness: bool, method: str = 'abs')[source]

Sets dipole corrector magnets to different setpoints

Sets horizontal or vertical CMs as specified in ords and skewness, respectively, to setpoints [rad] and updates the magnetic fields. If the corresponding setpoint exceeds the CM limit specified in the corresponding lattice field CMlimit, the CM is clipped to that value and a warning is being printed (to switch off, use warning(‘off’,’SC:CM1’)). Positive setpoints will result in kicks in the positive horizontal or vertical direction.

Parameters:
  • ords -- Array of CM ordinates in the lattice structure (ex: SC.ORD.CM[0])

  • setpoints -- CM setpoints (array or single value for all CMs) [rad]

  • skewness -- boolean specifying CM dimension ([False|True] -> [hor|ver])

  • method -- ‘abs’ (default), Use absolute setpoint ‘rel’, Use relative setpoint to current value ‘add’, Add setpoints to current value

Examples

Set all registered horizontal CMs to zero:

SC.set_cm_setpoints(ords=SC.ORD.HCM, skewness=False, setpoints=0.0)

Add 10urad to the fourth registered vertical CM:

SC.set_cm_setpoints(ords=SC.ORD.VCM[4], setpoints=1E-5, skewness=True, method='add')
set_magnet_setpoints(ords: int | List[int] | ndarray, setpoints: float | List[float] | ndarray, skewness: bool, order: int, method: str = 'abs', dipole_compensation: bool = False)[source]

Sets magnets to setpoints

Sets magnets (except CMs) as specified in ords to setpoints while order and skewness defines which field entry should be used (see below). The setpoints may be given relative to their nominal value or in absolute terms. If the considered quadrupole is a combined function magnet with non-zero bending angle and the kick compensation flag ‘dipole_compensation’=True, the appropriate bending angle difference is calculated and the horizontal CM setpoint is changed accordingly to compensate for that dipole kick difference. If the setpoint of a skew quadrupole exceeds the limit specified in the corresponding lattice field SkewQuadLimit, the setpoint is clipped to that value and a warning is logged.

Parameters:
  • ords -- Array of magnets ordinates in the lattice structure (ex: SC.ORD.HCM) (numpy.array() or list of int [int,int,..])

  • setpoints -- magnets setpoints (array or single value for all magnets). setpoints are assigned to the given order and skewness, i.e. once updated through SimulatedCommissioning.apply_errors, they correspond to a single element of PolynomA or PolynomB

  • skewness -- boolean specifying magnet plane ([False|True] -> [PolynomB|PolynomA])

  • method -- ‘abs’ (default), Use absolute setpoint ‘rel’, Use relative setpoint to nominal value ‘add’, Add setpoints to current value

  • order -- Numeric value defining the order of the considered magnet: [0,1,2,…] => [dip,quad,sext,…]

  • dipole_compensation -- (default = False) Used for combined function magnets. If this flag is set and if there is a horizontal CM registered in the considered magnet, the CM is used to compensate the bending angle difference if the applied quadrupole setpoints differs from the design value.

Examples

Identify the ordinates of all elements named ‘SF’ and switch their sextupole component off:

ords = sc_tools.ords_from_regex(SC.RING,'SF')
SC.register_magnets(ords)
SC.set_magnet_setpoints(ords=ords, skewness=False, order=2, setpoints=0.0, method='abs')

Identify the ordinates of all elements named QF and QD and set their quadrupole component to 99% of their design value:

ords = sc_tools.ords_from_regex(SC.RING,'QF|QD')
SC.register_magnets(ords)
SC.set_magnet_setpoints(ords=ords, skewness=False, order=1, setpoints=0.99, method='rel')
set_random_multipole_errors(ords: ndarray, BA)[source]

Applies random multipole errors specified in BA in the lattice elements ords of RING. It generates random multipole components with a 2-sigma truncated Gaussian distribution from each of the BA entries. The final multipole errors are stored in the PolynomA/BOffset of the lattice elements.

Parameters:
  • ords -- Ordinates of the considered magnets.

  • BA -- [N x 2] array of PolynomB/A multipole errors. [normal, skew] components.

Examples

Defines random multipole components for the ‘QF’ magnet and adds it to the field offsets of all magnets named ‘QF’:

ords = sc_tools.ords_from_regex(SC.RING,'QF')
BA = np.array([[1E-5, 0], [1E-4, 0], [0, 0], [1E-2, 0]])
SC.set_random_multipole_errors(ords, BA)

See also

SC.update_magnets, SC.set_systematic_multipole_errors

set_systematic_multipole_errors(ords: ndarray, BA, order: int, skewness: bool)[source]

Applies multipole errors specified in AB in the lattice elements ords of RING depending on the specified options. It sets the systematic multipoles of the field component defined by option ‘order’ and ‘type’. It is required that the BA entries are normalized by that component, e.g. BA[1, 0]=1 for skew-quadrupole systematic multipoles. The systematic multipoles are from now on scaled with the current magnet excitation and added to the PolynomA/B fields.

Parameters:
  • ords -- Ordinates of the considered magnets.

  • BA -- [N x 2] array of PolynomA/B multipole errors.

  • order -- Numeric value defining the order of the considered magnet: [0,1,2,…] => [dip,quad,sext,…]

  • skewness -- if True apply errors to skew fields (PolynomA) if False to normal fields (PolynomB)

Examples

Defines systematic multipole components for the ‘QF’ magnet and adds it to the field offsets of all magnets named ‘QF’:

ords = sc_tools.ords_from_regex(SC.RING,'QF')
BA = np.array([[1E-5, 0], [1E-4, 0], [0, 0], [1E-2, 0]])
RING = SC.set_systematic_multipole_errors(RING, ords, BA, 1, False)

See also

SC.update_magnets, SC.set_random_multipole_errors

support_offset_and_roll(s_locations: ndarray) Tuple[ndarray, ndarray][source]

This function evaluates the total offsets, roll, pitch and yaw angles of the support structures that have been defined via SC.register_support at the longitudinal positions s by linearly interpolating between support structure start- and endpoints (girder + sections + plinths, if registered). Note that this calculation may not provide the proper values if magnets with non-zero bending angle are within the support structure because it does not account for the rotation of the local coordinate system along the beam trajectory.

Parameters:

s_locations -- Array of s-positions at which the offset is evaluated.

Return type:

[3,length(s)]-array containing the [dx/dy/dz] total support structure offsets at s.

[3,length(s)]-array containing the [az/ax/ay] total support structure rolls at s.

See also

SC.register_support, SC.update_support, plot_support

update_cavities(ords: ndarray | None = None)[source]

Updates the cavity fields Voltage, Frequency and TimeLag in SC.RING as specified in ords. If no ordinates are given explicitly, all registered cavities defined in SC.ORD.RF are updated. For each cavity and each field, the setpoints, calibration errors and offsets are considered.

Parameters:

ords -- Cavity ordinates to be updated.

See also

SC.register_cavities, SC.apply_errors

update_magnets(ords: ndarray | None = None)[source]

Updates the magnets specified in RING as specified in ords. If no ordinates are given explicitly, all registered magnets defined in SC.ORD.Magnet are updated. For each magnet the setpoints (SetPointA/B) and calibration errors (CalErrorA/B) are evaluated. If systematic multipole components are specified, e.g. in SysPolBFromB for systematic PolynomB-multipoles induced by PolynomB entries, the corresponding multipole components are scaled by the current magnet excitation and added, as well as static field offsets (if specified in PolynomA/BOffset). If the considered magnet has a bending angle error (from pure bending angle error or due to a combined function magnet), the corresponding horizontal dipole magnetic field is calculated and added to the PolynomB(1) term. It is thereby assured that a dipole error doesn’t alter the coordinate system. If the considered magnet is registered as a slpit magnet (‘MasterOf’), the errors and setpoints of the master magnet are applied to the fields of the child magnets. Note that split quadrupole magnets with different gradients, however, or split CMs can currently not be updated correctly.

Parameters:

ords -- Magnets ordinates to be updated.

See also

SC.register_magnets, SC.apply_errors, SC.set_systematic_multipole_errors, SC.set_random_multipole_errors, set_magnet_setpoints, set_cm_setpoints

update_supports(offset_bpms: bool = True, offset_magnets: bool = True)[source]

This function updates the offsets and rolls of the elements in SC.RING based on the current support errors, by setting the lattice fields T1, T2, and R1, R2 for magnets and the fields SupportOffset and SupportRoll for BPMs.

Parameters:
  • offset_bpms -- If true, BPM offsets are updated.

  • offset_magnets -- If true, magnet offsets are updated.

See also

SC.register_support, SC.support_offset_and_roll, plot_support

verify_structure()[source]

Verifies the integrity of SimilatedCommissionining class instance and warns if things look fishy. If you find something that is missing please contact us.

See also

SC.register_magnets, SC.register_bpms, SC.register_cavities

Subclasses

This module contains the classes needed in the main data structure (SimulatedCommissioning).

class pySC.core.classes.DotDict(*args, **kwargs)[source]
deepcopy() DotDict[source]

Returns a deep copy

class pySC.core.classes.Indices[source]

Indices of magnets, RF cavities, BPMs, correctors, girders, plinths and sections in the AT lattice

Parameters:
  • BPM -- numpy.array of integers. index of beam position monitors

  • RF -- numpy.array of integers. index of RF cavities

  • Magnet -- numpy.array of integers. index of magnets

  • SkewQuad -- numpy.array of integers. index of skew quadrupoles

  • HCM -- numpy.array of integers. index of horizontal correctors

  • VCM -- numpy.array of integers. index of vertical correctors

  • Girder -- numpy.array of 2 element tuples. (start, end) start and end index of Girders

  • Plinth -- numpy.array of 2 element tuples. (start, end) start and end index of Plinths

  • Section -- numpy.array of 2 element tuples. (start, end) start and end index of Sections

Examples

get BPM indices:

SC.ORD.BPM
class pySC.core.classes.Injection[source]

Define injection parameters for “pySC”

properties of this class are:
beamLostAt:

(default=1) Relative amount of particles which may be lost before BPM reading is NaN

Z0ideal:

(default=numpy.zeros(6)) Design injected trajectory

Z0:

(default=numpy.zeros(6)) Injected trajectory

beamSize:

(default=numpy.zeros((6,6))) Injected bunch beam size

randomInjectionZ:

(default=numpy.zeros(6)) Injected beam random trajectory jitter

staticInjectionZ:

(default=numpy.zeros(6)) Injected beam static trajectory offset

nParticles:

(default=1) Number of particles per bunch

nTurns:

(default=1) Number of turns for tracking

nShots:

(default=1) Number of injections for averaging BPM reading

trackMode:
(default=’TBT’) Tracking mode can be one of:

TBT (Turn By Turn), ORB (Closed Orbit), PORB (pseudo Closed Orbit) turn-by-turn tracking with trajectories averaged over the turns

class pySC.core.classes.Sigmas[source]

Defined rms (sigma) size of errors later assigned to BPM, RF, magnets, supports and injection

Parameters:
  • BPM -- errors sigma’s in BPMs

  • Magnet -- errors sigma’s in Magnets

  • RF -- errors sigma’s in RF

  • Support -- errors sigma’s in Supports

  • randomInjectionZ -- numpy.array of 6 elements defining sigmas of random errors for the input coordinates x, x’, y, y’, delta and ct

  • staticInjectionZ -- numpy.array of 6 elements defining an offset for the input coordinates x, x’, y, y’, delta and ct

  • Circumference -- float. Circumference error

Examples

get BPM sigmas from SimualtedCommissioning:

SC.SIG.BPM