Correction algorithms

Trajectory and Orbit

This module contains functions to correct trajectory (first turns(s)) and orbit. In all functions, the provided response matrix needs to be matching:

registered BPMs in SC.ORD.BPM or bpm_ords if provided registered CMs in SC.ORD.CM or cm_ords if provided assumes all BPMs are used both horizontal and vertical plane

pySC.correction.orbit_trajectory.balance(SC, response_matrix, reference=None, cm_ords=None, bpm_ords=None, eps=0.0001, maxsteps=10, **pinv_params)[source]

Balance two-turn BPM readings

Generates a period-1 closed orbit, after two-turn transmission has been achieved. This is done by iteratively applying correction steps, calculated based on the pseudo-inverse two-turn trajectory response matrix Mplus. The trajectory in the first turn is corrected towards the reference orbit reference, whereas the trajectory in the second turn is corrected towards the trajectory measured in the first turn; this approach seems to be more stable than the directly application of the two-turn TRM to the two-turn BPM readings. It converges to a state where BPM readings in both turns are very similar, indicating a period-1 closed orbit.

Parameters:
  • SC -- SimulatedCommissioning class instance.

  • response_matrix -- Trajectory-response matrix.

  • reference -- (None) target orbit in the format [x_1 … x_n y_1 …y_n], where [x_i,y_i]` is the target position at the i-th BPM.

  • cm_ords -- List of CM ordinates to be used for correction (SC.ORD.CM)

  • bpm_ords -- List of BPM ordinates at which the reading should be evaluated (SC.ORD.BPM)

  • maxsteps -- break, if this number of correction steps have been performed (default = 100)

  • eps -- break, if the coefficient of variation of the RMS BPM reading is below this value

Returns:

SimulatedCommissioning class instance with corrected SC.RING

pySC.correction.orbit_trajectory.correct(SC, response_matrix, reference=None, cm_ords=None, bpm_ords=None, eps=0.0001, target=0, maxsteps=30, scaleDisp=0, **pinv_params)[source]

Iterative orbit/trajectory correction

Iteratively applies orbit corrections using the pseudoinverse of the trajectory response_matrix, until a break-condition specified by one of ‘eps’, ‘target’ or ‘maxsteps’ is met. The dispersion can be included, thus the rf frequency as a correction parameter. If the dispersion is to be included, response_matrix has to have the size (2 * len(SC.ORD.BPM)) x (len(SC.ORD.HCM) + len(SC.ORD.VCM) + 1), otherwise the size (2 * len(SC.ORD.BPM)) x (len(SC.ORD.HCM) + len(SC.ORD.VCM)), or correspondingly if the CM and/or BPM ordinates for the correction is explicitly given (see options below). SC.RING is assumed to be a lattice with transmission through all considered turns. Raises RuntimeError if transmission is lost.

Parameters:
  • SC -- SimulatedCommissioning class instance.

  • response_matrix -- Trajectory/orbit-response matrix.

  • reference -- (None) target orbit in the format [x_1 … x_n y_1 …y_n], where [x_i,y_i]` is the target position at the i-th BPM.

  • cm_ords -- List of CM ordinates to be used for correction (SC.ORD.CM)

  • bpm_ords -- List of BPM ordinates at which the reading should be evaluated (SC.ORD.BPM)

  • maxsteps -- break, if this number of correction steps have been performed (default = 100)

  • eps -- break, if the coefficient of variation of the RMS BPM reading is below this value

  • target -- (default =0 ) break, if the RMS BPM reading reaches this value

  • scaleDisp -- (default =0 ) Scaling factor for and flag indicating if the dispersion is included in the response matrix

Returns:

SimulatedCommissioning class instance with corrected SC.RING

Examples

Switch to orbit mode, get the model response matrix and dispersion. Calculate the psudo-inverse while scaling the dispersion by 1E7 and using a Tikhonov regularization parameter of 10. Finally, apply and apply orbit feedback including dispersion:

SC.INJ.trackMode = 'ORB'
MCO = SCgetModelRM(SC, SC.ORD.BPM, SC.ORD.CM, trackMode='ORB')
eta = SCgetModelDispersion(SC, SC.ORD.BPM, SC.ORD.RF,
                           trackMode='ORB', Z0=np.zeros(6),
                           nTurns=1, rfStep=1E3,
                           useIdealRing=True)
SC = correct(SC, np.column_stack((MCO, 1E8 * eta)), alpha=10, target=0, maxsteps=50, scaleDisp=1E8)
pySC.correction.orbit_trajectory.first_turn(SC, response_matrix, reference=None, cm_ords=None, bpm_ords=None, maxsteps=100, **pinv_params)[source]

Achieves one-turn transmission

Achieves a first turn in SC.RING. This algorithm is based on the idea that repeated trajectory corrections calculated via a suitably regularized pseudo-inverse of response matrix will drive the BPM readings and CM settings to a fixed point.

lim_{n->oo} B_n = const. , with B_{n+1} = Phi(response_matrix^{-1} . B_{n} ), (1)

where mapping Phi maps corrector kicks to BPM-readings. The RMS-values of both, BPM readings and CM settings, are determined by the regularization of pseudo-inverted response matrix. Successively - during the course of repeated application of (1) - more and more transmission is achieved throughout the ring, more magnets are traversed near their magnetic center (which is hopefully at least somewhere near the BPM zero-point), resulting in decreased kicks. Otherwise,if trajectory correction cannot proceed further to next BPM the kicks of an increasing number of the last reached CMs are deterministically ``wiggled’’ until transmission to the next BPM is achieved. Then, application of (1) is resumed.

Parameters:
  • SC -- SimulatedCommissioning class instance.

  • response_matrix -- Trajectory-response matrix.

  • reference -- (None) target orbit in the format [x_1 … x_n y_1 …y_n], where [x_i,y_i]` is the target position at the i-th BPM.

  • cm_ords -- List of CM ordinates to be used for correction (SC.ORD.CM)

  • bpm_ords -- List of BPM ordinates at which the reading should be evaluated (SC.ORD.BPM)

  • maxsteps -- break, if this number of correction steps have been performed (default = 100)

Returns:

SimulatedCommissioning class instance with corrected SC.RING

pySC.correction.orbit_trajectory.stitch(SC, response_matrix, reference=None, cm_ords=None, bpm_ords=None, n_bpms=4, maxsteps=30, **pinv_params)[source]

Achieves 2-turn transmission

The purpose of this function is to go from 1-turn transmission to 2-turn transmission. This is done by applying orbit correction based on the pseudo inverse trajectory response matrix ‘Mplus’ applied to the first BPMs in the ‘SC.RING’. The reading of the BPMs in the second turn is corrected towards the reading of these BPMs in the first turn. This approach has been seen to be more stable than the direct application of the two-turn inverse response matrix to the two-turn BPM data.

Parameters:
  • SC -- SimulatedCommissioning class instance.

  • response_matrix -- Trajectory-response matrix.

  • reference -- (None) target orbit in the format [x_1 … x_n y_1 …y_n], where [x_i,y_i]` is the target position at the i-th BPM.

  • cm_ords -- List of CM ordinates to be used for correction (SC.ORD.CM)

  • bpm_ords -- List of BPM ordinates at which the reading should be evaluated (SC.ORD.BPM)

  • n_bpms -- Number of BPMs to match the trajectory in first and second turn

  • maxsteps -- break, if this number of correction steps have been performed (default = 100)

Returns:

SimulatedCommissioning class instance with corrected SC.RING

Examples

Calculate the 2-turn response matrix and get the pseudo inverse using a Tikhonov regularization parameter of 10. Switch the injection pattern to 2 turns and apply the stitching using the first three BPMs, for a maximum of 20 steps:

RM2 = SCgetModelRM(SC, SC.ORD.BPM, SC.ORD.CM, nTurns=2)
SC.INJ.nTurns = 2
SC = stitch(SC, RM2, alpha=10, nBPMs=3, maxsteps=20)

RF

This module contains functions to correct RF phase and frequency.

pySC.correction.rf.phase_and_energy_error(SC, cav_ords)[source]

This is not a simple observable in reality.

Tune

This module contains functions to correct betatron tunes.

pySC.correction.tune.tune_scan(SC, quad_ords, rel_quad_changes, target=1, n_points=60, do_plot=False, nParticles=None, nTurns=None, full_scan=False)[source]

Varies two quadrupole families to improve beam transmission, on a grid of relative setpoints specified in rel_quad_changes in a spiral-like pattern to increase the beam transmission.

Parameters:
  • SC -- SimulatedCommissioning instance

  • quad_ords -- [1x2] array of quadrupole ordinates {[1 x NQ1],[1 x NQ2]}

  • rel_quad_changes -- [1x2] cell array of quadrupole setpoints {[SP1_1,…,SP1_N1],[SP2_1,…,SP2_N2]} with N2=N1

  • target (int, optional) -- Transmission target at nTurns

  • n_points (int, optional) -- Number of points for the scan

  • nParticles (int, optional) -- Number of particles used for tracking (for convenience, otherwise SC.INJ.nParticles is used)

  • nTurns (int, optional) -- Number of turns used for tracking (for convenience, otherwise SC.INJ.nTurns is used)

  • full_scan (bool , optional) -- If false, the scan finishes as soon as the target is reached

  • do_plot (bool , optional) -- If true, beam transmission is plotted at every step

Returns:

Updated SC structure with applied setpoints for maximised beam transmission Relative setpoints which satisfied the target condition if reached, or the values which resulted in best transmission Number of achieved turns Turn-by-turn particle loss

see also: bpm_reading, generate_bunches

Chromaticity

This module contains functions to fit the chromaticity of ‘SC.RING’.

pySC.correction.chroma.fit_chroma(SC, s_ords, target_chroma=None, init_step_size=array([2, 2]), xtol=0.0001, ftol=0.001)[source]

Applies a chromaticity correction using two sextupole families.

Parameters:
  • SC -- SimulatedCommissioning instance

  • s_ords -- [2xN] array or list [[1 x NS1],[1 x NS2], [1 x NS3], …] of sextupole ordinates

  • target_chroma ([1x2] array, optional) -- Target chromaticity for correction. Default: chromaticity of ‘SC.IDEALRING’

  • init_step_size ([1x2] array, optional) -- Initial step size for the solver. Default: [2,2]

  • xtol (float, optional) -- Step tolerance for solver. Default: 1e-4

  • ftol (float, optional) -- Merit tolerance for solver. Default: 1e-3

Returns:

SimulatedCommissioning instance with corrected chromaticity.

Return type:

SC

Example

SC = fit_chroma(SC, s_ords=[SCgetOrds(sc.RING, ‘SF’), SCgetOrds(sc.RING, ‘SD’)], target_chroma=numpy.array([1,1]))