Constants

These constants rely a lot on the HL-LHC/LHC Naming Scheme.

Equation System

Builds and solves the equation system from the rdt-to-corrector-maps given. This is Eq. (40) of [DillyNonlinearIRCorrections2023].

irnl_rdt_correction.equation_system.build_equation_system(rdt_maps: Sequence[Dict[RDT, Sequence[str]]], correctors: Sequence[IRCorrector], ip: int, optics_seq: Sequence[Optics], feeddown: int) Tuple[Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes], Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes]][source]

Builds equation system as in Eq. (40) of [DillyNonlinearIRCorrections2023] for a given ip for all given optics and error files (i.e. beams) and grouped RDTs, i.e. RDTs that share correctors.

Parameters:
  • rdt_maps (Sequence[RDTMap]) -- Sequence of RDTMap objects, i.e. a dict mapping the RDT to the used correctors. This is a subset of all RDTs, as they have been grouped by common correctors before.

  • correctors (Sequence[IRCorrector]) -- IRCorrectors to be used.

  • ip (int) -- Current IP to correct.

  • optics_seq (Sequence[Optics]) -- Sequence of given Optics (twiss and errors)

  • feeddown (int) -- Orders of feed-down to include calculating the integral on the rhs

Returns:

tuple of

b_matrix: np.array N_rdts x N_correctors integral: np.array N_rdts x 1

irnl_rdt_correction.equation_system.get_available_correctors(field_components: Set[str], accel: str, ip: int, optics_seq: Sequence[Optics]) List[IRCorrector][source]

Gets the available correctors by checking for their presence in the optics. If the corrector is not found in this ip, the ip is skipped. If only one corrector (left or right) is present a warning is issued. If one corrector is present in only one optics (and not in the other) an Environment Error is raised.

Parameters:
  • field_components (Set[str]) -- Set of field components to be corrected (i.e. correctors to be found for). Field components are e.g. “a3”, “b4” etc.

  • accel (str) -- Which accelerator to use.

  • ip (int) -- Which IP we are currently working on.

  • optics_seq (Sequence[Optics]) -- Sequence of Optics (twiss and errors).

Returns:

List of IRCorrectors to use to correct given field components.

Return type:

List[IRCorrector]

irnl_rdt_correction.equation_system.get_corrector_coefficient(rdt: RDT, corrector: IRCorrector, optics: Optics) float[source]

Calculate B-Matrix Element for Corrector. These are the entries on the lhs of Eq. (40) of [DillyNonlinearIRCorrections2023] including feed-down coefficient z and signs. Any imaginary i coefficient is also included, making all values real.

Parameters:
  • rdt (RDT) -- Current RDT

  • corrector (IRCorrector) -- IRCorrector at which the coefficient is to be calculated.

  • optics (Optics) -- Current optics

Returns:

Calculated matrix value.

Return type:

float

irnl_rdt_correction.equation_system.get_current_rdt_maps(rdt_maps: Sequence[Dict[RDT, Sequence[str]]]) Tuple[Sequence[Dict[RDT, Sequence[str]]], Sequence[Dict[RDT, Sequence[str]]], Set[str]][source]

Creates a new rdt_map with all rdt’s that share correctors.

This function is called in a while-loop, so rdt_maps is the remaining_rdt_maps from the last loop. The while-loop is interrupted when no remaining rdts are left.

Parameters:

rdt_maps (Sequence[RDTMap]) -- (Still) available RDTMaps to be checked.

Returns:

3-Tuple consisting of the sequence of current RDTMaps to use, the remaining RDT maps and the checked (i.e. to be used) correctors in this loop.

Return type:

Tuple[Sequence[RDTMap], Sequence[RDTMap], Set[str]]

irnl_rdt_correction.equation_system.get_elements_integral(rdt: RDT, ip: int, optics: Optics, feeddown: int) float[source]

Calculate the RDT integral for all elements of the IP. These are the entries on the rhs of Eq. (40) of [DillyNonlinearIRCorrections2023], including sign.

Parameters:
  • rdt (RDT) -- Current RDT

  • ip (int) -- Current IP

  • optics (Optics) -- Current optics

  • feeddown (int) -- Orders of feed-down to include

Returns:

Calculated Integral value.

Return type:

float

irnl_rdt_correction.equation_system.get_side_sign(n: int, side: str) int[source]

Sign of the integral and corrector for this side.

This is the exp(iπnθ(s_w−s_IP)) part of Eq. (16) or (-1)^(θ(s_w−s_IP)) of later equations, e.g. Eq. (17) in [DillyNonlinearIRCorrections2023].

Parameters:
  • n (int) -- order of the RDT

  • side (str) -- side of the IP

Returns:

Either -1 or 1

Return type:

int

irnl_rdt_correction.equation_system.init_corrector_and_optics_values(correctors: Sequence[IRCorrector], optics_seq: Sequence[Optics], update_optics: bool, ignore_settings: bool) Dict[IRCorrector, Sequence[float]][source]

Save original corrector values from optics (for later restoration, only if update_optics is False) and if ignore_settings is True, the corrector values in the optics are set to 0. Otherwise, the corrector object value is initialized with the value from the optics. An error is thrown if the optics differ in value.

Parameters:
  • correctors (Sequence[IRCorrector]) -- Sequence of IRCorrectors to initialize.

  • optics_seq (Sequence[Optics]) -- Optics to get values from.

  • update_optics (bool) -- If not set, saves initial data for later recovery.

  • ignore_settings (bool) -- If set, values in Optics will be set to zero for given correctors.

Returns:

The saved values per corrector (per Optics).

These are the values directly taken from the optics, i.e. the sign might be different to the sign stored in the correctors objects. If the optics are updated anyway, this is an empty dict.

Return type:

Dict[IRCorrector, Sequence[float]]

irnl_rdt_correction.equation_system.optics_update(correctors: Sequence[IRCorrector], optics_seq: Sequence[Optics]) None[source]

Updates the corrector strength values in the current optics.

Parameters:
  • correctors (Sequence[IRCorrector]) -- IRCorrectors containing the new values.

  • optics_seq (Sequence[Optics]) -- Optics to update.

irnl_rdt_correction.equation_system.restore_optics_values(saved_values: Dict[str, Sequence[float]], optics_seq: Sequence[Optics]) None[source]

Restore saved initial corrector values (if any) to optics.

Parameters:
  • saved_values (Dict) -- Saved initial values (key: entry in optics, values: value per optics)

  • optics_seq (Sequence[Optics]) -- Optics to update.

irnl_rdt_correction.equation_system.solve(rdt_maps: Sequence[Dict[RDT, Sequence[str]]], optics_seq: Sequence[Optics], accel: str, ips: Sequence[int], update_optics: bool, ignore_corrector_settings: bool, feeddown: int, iterations: int, solver: str) Sequence[IRCorrector][source]

Calculate corrections, i.e. build and solve Eq. (40) of [DillyNonlinearIRCorrections2023]. Corrections are performed by grouping RDTs with common correctors. If possible, these are ordered from the highest order to lowest, to be able to update optics and include their feed-down.

Parameters:
  • rdt_maps (Sequence[RDTMap]) -- Sequence of RDTMap objects, i.e. a dict mapping the RDT to the used correctors.

  • optics_seq (Sequence[Optics]) -- Sequence of Optics objects containing twiss and errors.

  • accel (str) -- Accelerator to use (implemented ‘lhc’ and ‘hllhc’).

  • ips (Sequence[int]) -- Sequence of IPs to correct. Elements will split by IP, assuming their name ends with “L” or “R” for left and right followed by the IP number.

  • update_optics (bool) -- Sorts the RDTs to start with the highest order and updates the corrector strengths in the optics after calculation, so the feeddown to lower order correctors is included.

  • ignore_corrector_settings (bool) -- Ignore the current settings of the correctors. If this is not set the corrector values of the optics are used as initial conditions.

  • feeddown (int) -- Orders of feed-down to include calculating the integral on the rhs of Eq. (40) of [DillyNonlinearIRCorrections2023]

  • iterations (int) -- (Re-)iterate correction, starting with the previously calculated values. Needs to be > 0, as the first calculation counts as an iteration.

  • solver (str) -- Solver to use: ‘inv’, ‘lstq’ or ‘linear’

Returns:

Sequence of IRCorrector objects, which define the

IR corrector and its value and contain additional information about the corrector.

Return type:

Sequence[IRCorrector]

irnl_rdt_correction.equation_system.solve_equation_system(correctors: Sequence[IRCorrector], lhs: array, rhs: array, solver: str)[source]

Solves the system with the given solver.

The result is transferred to the corrector values internally.

Parameters:
  • correctors (Sequence[IRCorrector]) -- Sequence of IRCorrectors this equation system was build for. Stores the resulting values.

  • lhs (np.array) -- Left hand side of the Eqs to solve.

  • rhs (np.array) -- Right hand side of the Eqs to solve.

  • solver (str) -- Solver to use (see SOLVER_MAP below).

Input Options

Handles the input parameters, contains their default values and checks for their validity.

class irnl_rdt_correction.input_options.InputOptions(beams: Sequence[int], twiss: Sequence[str | Path | DataFrame | TfsDataFrame], errors: Sequence[str | Path | DataFrame | TfsDataFrame | None] = None, rdts: Sequence[str] = None, rdts2: Sequence[str] = None, accel: str = 'lhc', feeddown: int = 0, ips: Sequence[int] = (1, 2, 5, 8), solver: str = 'lstsq', update_optics: bool = True, iterations: int = 1, ignore_corrector_settings: bool = False, ignore_missing_columns: bool = False, output: str = None)[source]

DataClass to store the input options. On creation it asserts that the input parameters make sense and adds what’s missing. Checks include:

  • Check accelerator name is valid

  • Set default RDTs if not given (see DEFAULT_RDTS)

  • Check required parameters are present (twiss, errors, beams, rdts)

  • Check feeddown and iterations

classmethod from_args_or_dict(opt: dict | InputOptions | None = None) InputOptions[source]

Create an InputOptions instance from the given dictionary. If the input is empty, arguments will be parsed from commandline.

Parameters:

opt (Union[dict, DotDict]) -- Function options in dictionary format. Description of the arguments are given in irnl_rdt_correction.main.irnl_rdt_correction(). Optional, if not given parses commandline args

Returns:

(Parsed and) checked options.

Return type:

InputOptions

irnl_rdt_correction.input_options.allow_commandline_and_kwargs(func)[source]

Decorator to allow a function to take options from the commandline or via kwargs, or given an InputOptions instance.

IO Handling

Functions for reading input and writing output.

irnl_rdt_correction.io_handling.build_correction_df(correctors: Sequence[IRCorrector]) TfsDataFrame[source]

Build the DataFrame that contains the information for corrector powering.

Parameters:

correctors (Sequence[IRCorrector]) -- Sequence or IRCorrectors containing their desired settings.

Returns:

Converted DataFrame

Return type:

TfsDataFrame

irnl_rdt_correction.io_handling.build_correction_str(correctors: Sequence[IRCorrector]) str[source]

Build the MAD-X command that contains the information for corrector powering.

Parameters:

correctors (Sequence[IRCorrector]) -- Sequence or IRCorrectors containing their desired settings.

Returns:

Converted MAD-X commands.

Return type:

str

irnl_rdt_correction.io_handling.build_correction_str_from_df(correctors_df: DataFrame) str[source]

Build the MAD-X command that contains the information for corrector powering from the given DataFrame.

Parameters:

correctors_df (DataFrame) -- Corrector information in DataFrame format.

Returns:

Converted MAD-X commands.

Return type:

str

irnl_rdt_correction.io_handling.convert_numeric_columns_to_float(df: TfsDataFrame) TfsDataFrame[source]

Convert numeric columns in df to float. This avoids that the user accidentally gives columns with dtype int (e.g. all 0), in which case assigning float values will fail in some pandas version after 2.1.0 (which had a deprecation warning).

irnl_rdt_correction.io_handling.get_and_write_output(out_path: str | Path, correctors: Sequence[IRCorrector]) Tuple[str, TfsDataFrame][source]

Convert the IRCorrector data to MAD-X commands and TFSDataFrame and write to files.

Parameters:
  • out_path (Union[str, Path]) -- Path to output file (suffix will be added/changed to fit format)

  • correctors (Sequence[IRCorrector]) -- Sequence or IRCorrectors containing their desired settings.

Returns:

IRCorrectors information converted to MAD-X command and TFSDataFrame

Return type:

Tuple[str, tfs.TfsDataFrame]

irnl_rdt_correction.io_handling.get_optics(beams: Sequence[int], twiss: Sequence[str | Path | DataFrame | TfsDataFrame], errors: Sequence[str | Path | DataFrame | TfsDataFrame | None], orders: Sequence[int], ignore_missing_columns: bool) Tuple[Optics][source]

Get the Optics instances from beams, twiss and errors. Also checks the DataFrames for containing needed information.

Parameters:
  • beams (Sequence[int]) -- Beam number the twiss/errors refer to

  • twiss (Sequence[StrOrPathOrDataFrame]) -- Sequence of twiss TfsDataFrames or paths to the tfs files.

  • errors (Sequence[StrOrPathOrDataFrame]) -- Sequence of errors TfsDataFrames or paths to the tfs files. Can be None.

  • orders (Sequence[int]) -- Orders needed for calculations (used to check DataFrames)

  • ignore_missing_columns (bool) -- If set, fills missing columns with zeros.

Returns:

Tuple of Optics objects containing combined information

about beam, twiss and errors.

Return type:

Sequence[Optics]

irnl_rdt_correction.io_handling.get_tfs(paths: Sequence[str | Path | DataFrame | TfsDataFrame]) Sequence[DataFrame][source]

Get TFS-DataFrames from Paths or return DataFrames if already loaded.

Parameters:

paths (Sequence[StrOrPathOrDataFrame]) -- (TFS)DataFrames or paths to tfs files.

Returns:

Loaded (TFS)DataFrames.

Return type:

Sequence[DataFrame]

irnl_rdt_correction.io_handling.maybe_switch_signs(optics: Optics)[source]

Switch the signs for Beam optics. This is due to the switch in direction between beam and (anti-) symmetry after a rotation of 180deg around the y-axis of magnets. Magnet orders that show anti-symmetry are: a1 (K0SL), b2 (K1L), a3 (K2SL), b4 (K3L) etc.

This brings the Beam 2 KNL and x values to Beam 4 definition.

RDT Handling

Class definitions of RDTs and IRCorrectors and functions handling these classes for checking and sorting.

class irnl_rdt_correction.rdt_handling.IRCorrector(field_component: str, accel: str, ip: int, side: str)[source]

Representation of an IR-Corrector in the accelerator.

class irnl_rdt_correction.rdt_handling.RDT(name: str)[source]

Representation of Resonance Driving Term attributes.

irnl_rdt_correction.rdt_handling.check_corrector_order(rdt_maps: Sequence[Dict[RDT, Sequence[str]]], update_optics: bool, feed_down: int)[source]

Perform checks on corrector orders compared to RDT orders and feed-down.

Parameters:
  • rdt_maps (Sequence[RDTMap]) -- RDTMaps to check.

  • update_optics (bool) -- True if optics should be updated after each iteration.

  • feed_down (int) -- Order of feed-down to include.

irnl_rdt_correction.rdt_handling.get_needed_orders(rdt_maps: Sequence[Dict[RDT, Sequence[str]]], feed_down: int) Sequence[int][source]

Returns the sorted orders needed for correction, based on the order of the RDTs to correct plus the feed-down involved and the order of the corrector, which can be higher than the RDTs in case one wants to correct via feed-down.

Parameters:
  • rdt_maps (Sequence[RDTMap]) -- RDTMaps to check

  • feed_down (int) -- Order of maximum feed-down to include

Returns:

Needed field orders to calculate all corrections including feed-down

(sorted low-to-high).

Return type:

Sequence[int]

irnl_rdt_correction.rdt_handling.sort_rdts(rdts: Sequence[str] | Dict[str, Sequence[str]], rdts2: Sequence[str] | Dict[str, Sequence[str]]) Tuple[Dict[RDT, Sequence[str]], Dict[RDT, Sequence[str]]][source]

Sorts RDTs by reversed-order (i.e. high-to-low) and orientation (skew, normal).

Parameters:
  • rdts (RDTInputTypes) -- RDTs for first optics

  • rdts2 (RDTInputTypes) -- RDTs for second optics. If not given first optics RDTs are assumed.

Returns:

RDT mapping (RDT to corrector fields) for both optics.

Return type:

Tuple[RDTMap, RDTMap]

Utilities

Additional utilities used throughout the correction.

class irnl_rdt_correction.utilities.Optics(beam: int, twiss: DataFrame, errors: DataFrame)[source]

Store Optics Data.

class irnl_rdt_correction.utilities.Timer(name: str = 'start', print_fun: ~typing.Callable[[str], None] = <built-in function print>)[source]

Collect Times and print a summary at the end.

irnl_rdt_correction.utilities.corrector_sign_beam_symmetry(beam: int, column_name: str) int[source]

This function returns -1 if there is if there is a sign change change for the powering of the corrector and it’s field component, due to the anti-symmetric magnetic field. See Chapter 3.1 Beam Directions in [DillyNonlinearIRCorrections2023]

Parameters:
  • beam (int) -- beam used

  • column_name (str) -- column name of the field strength

Returns:

-1 if the corrector circuit has opposite sign to its field, 1 otherwise.

Return type:

int

irnl_rdt_correction.utilities.get_max_knl_order(df: TfsDataFrame) int[source]

Return the maximum order in the KN(S)L columns of the DataFrame.

irnl_rdt_correction.utilities.i_pow(n: int) complex[source]

i to the power of n.

irnl_rdt_correction.utilities.is_anti_mirror_symmetric(column_name: str) bool[source]

Returns true if the column name is a KNL/KNSL column and the magnetic field that this column represents is anti-symmetric upon mirroring on y axis.

irnl_rdt_correction.utilities.log_setup()[source]

Set up a basic logger.