Skip to content

Commit 8c7cd3a

Browse files
authored
Fixes main import and deprecation warnings (#14)
* fix absolute import * update gitignore * Silence rcond warning in lstsq in test * use math instead of np.math * added commandline tests * bumpversion
1 parent f8104ed commit 8c7cd3a

File tree

9 files changed

+177
-19
lines changed

9 files changed

+177
-19
lines changed

.gitignore

Lines changed: 4 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -260,4 +260,7 @@ Temporary Items
260260
*.out
261261
*.synctex.*
262262
*.toc
263-
*.fdb_latexmk
263+
*.fdb_latexmk
264+
265+
# Personal testing scripts
266+
tst_*

irnl_rdt_correction/__init__.py

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -13,7 +13,7 @@
1313
__title__ = "irnl-rdt-correction"
1414
__description__ = "Correction script to power the nonlinear correctors in the (HL-)LHC insertion regions based on RDTs."
1515
__url__ = "https://github.com/pylhc/irnl_rdt_correction"
16-
__version__ = "1.1.1"
16+
__version__ = "1.1.2"
1717
__author__ = "pylhc"
1818
__author_email__ = "[email protected]"
1919
__license__ = "MIT"

irnl_rdt_correction/__main__.py

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1,5 +1,5 @@
11
from irnl_rdt_correction.main import irnl_rdt_correction
2-
from utilities import log_setup
2+
from irnl_rdt_correction.utilities import log_setup
33

44
if __name__ == '__main__':
55
log_setup()

irnl_rdt_correction/equation_system.py

Lines changed: 3 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -7,6 +7,7 @@
77
88
"""
99
import logging
10+
import math
1011
from typing import Dict, List, Sequence, Set, Tuple
1112

1213
import numpy as np
@@ -354,7 +355,7 @@ def get_elements_integral(rdt: RDT, ip: int, optics: Optics, feeddown: int) -> f
354355
iksl_err = 1j*errors_df.loc[elements, f"K{n_mad:d}SL"]
355356

356357
k_sum += ((kl_opt + kl_err + iksl_opt + iksl_err) *
357-
(dx_idy**q) / np.math.factorial(q))
358+
(dx_idy**q) / math.factorial(q))
358359

359360
# note the minus sign before the sum!
360361
integral += -sum(np.real(i_pow(lm) * k_sum.to_numpy()) * (side_sign * betax * betay).to_numpy())
@@ -406,7 +407,7 @@ def get_corrector_coefficient(rdt: RDT, corrector: IRCorrector, optics: Optics)
406407
dx = twiss_df.loc[corrector.name, X] + errors_df.loc[corrector.name, f"{DELTA}{X}"]
407408
dy = twiss_df.loc[corrector.name, Y] + errors_df.loc[corrector.name, f"{DELTA}{Y}"]
408409
dx_idy = dx + 1j*dy
409-
z_cmplx = (dx_idy**p) / np.math.factorial(p) # Eq. (32)
410+
z_cmplx = (dx_idy**p) / math.factorial(p) # Eq. (32)
410411

411412
# Get the correct part of z_cmplx, see Eq. (36) in [DillyNonlinearIRCorrections2023]_
412413
if (corrector.skew and is_odd(lm)) or (not corrector.skew and is_even(lm)):

irnl_rdt_correction/io_handling.py

Lines changed: 16 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -220,8 +220,9 @@ def _check_dfs(beams: Sequence[int], twiss_dfs: Sequence[DataFrame], errors_dfs:
220220
- All needed field strengths are present in twiss
221221
-> either Error or Warning depending on ``ignore_missing_columns``.
222222
"""
223-
twiss_dfs, errors_dfs = tuple(tuple(df.copy() for df in dfs)
224-
for dfs in (twiss_dfs, errors_dfs))
223+
twiss_dfs, errors_dfs = tuple(tuple(
224+
convert_numeric_columns_to_float(df) for df in dfs) for dfs in (twiss_dfs, errors_dfs)
225+
)
225226

226227
needed_twiss = list(PLANES)
227228
needed_errors = [f"{DELTA}{p}" for p in PLANES]
@@ -258,7 +259,7 @@ def _check_dfs(beams: Sequence[int], twiss_dfs: Sequence[DataFrame], errors_dfs:
258259
f"the errors: {list2str(not_found_optics.to_list())}."
259260
f"They are assumed to be zero.")
260261
for indx in not_found_optics:
261-
errors.loc[indx, :] = 0
262+
errors.loc[indx, :] = 0.0
262263

263264
for df, needed_cols, file_type in ((
264265
twiss, needed_twiss, "twiss"), (errors, needed_errors, "error")
@@ -273,10 +274,21 @@ def _check_dfs(beams: Sequence[int], twiss_dfs: Sequence[DataFrame], errors_dfs:
273274
raise IOError(text)
274275
LOG.warning(text + " They are assumed to be zero.")
275276
for kl in missing_cols:
276-
df[kl] = 0
277+
df[kl] = 0.0
277278
return beams, twiss_dfs, errors_dfs
278279

279280

281+
def convert_numeric_columns_to_float(df: TfsDataFrame) -> TfsDataFrame:
282+
""" Convert numeric columns in df to float.
283+
This avoids that the user accidentally gives columns with dtype int (e.g. all 0),
284+
in which case assigning float values will fail in some pandas version after 2.1.0
285+
(which had a deprecation warning). """
286+
df = df.copy()
287+
numeric_columns = df.select_dtypes(include=['number']).columns
288+
df[numeric_columns] = df[numeric_columns].astype('float64')
289+
return df
290+
291+
280292
def maybe_switch_signs(optics: Optics):
281293
""" Switch the signs for Beam optics.
282294
This is due to the switch in direction between beam and

tests/helpers.py

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -60,7 +60,7 @@ def read_lhc_model(beam: int) -> TfsDataFrame:
6060

6161

6262
def generate_pseudo_model(n_ips: int, n_magnets: int, accel: str,
63-
betax: float = 1, betay: float = 1, x: float = 0, y: float = 0) -> pd.DataFrame:
63+
betax: float = 1.0, betay: float = 1.0, x: float = 0.0, y: float = 0.0) -> pd.DataFrame:
6464
"""Generate a Twiss-Like DataFrame with magnets as index and Beta and Orbit columns."""
6565
df = pd.DataFrame(
6666
index=(
@@ -77,7 +77,7 @@ def generate_pseudo_model(n_ips: int, n_magnets: int, accel: str,
7777
return df
7878

7979

80-
def generate_errortable(index: Sequence, value: float = 0, max_order: int = MAX_N) -> pd.DataFrame:
80+
def generate_errortable(index: Sequence, value: float = 0.0, max_order: int = MAX_N) -> pd.DataFrame:
8181
"""Return DataFrame from index and KN(S)L + D[XY] columns."""
8282
return pd.DataFrame(value,
8383
index=index,

tests/test_commandline.py

Lines changed: 141 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,141 @@
1+
import subprocess
2+
import sys
3+
from pathlib import Path
4+
5+
import pytest
6+
import tfs
7+
8+
from tests.helpers import (CIRCUIT, EPS, IP, MAX_N, STRENGTH, VALUE, generate_errortable,
9+
generate_pseudo_model, get_some_magnet_names)
10+
11+
12+
def test_package_execution_no_arguments():
13+
""" Tests if the package can be run as a module.
14+
Without arguments it should abort with an error message stating required arguments. """
15+
# Run the package as a module using subprocess
16+
result = subprocess.run([sys.executable, "-m", "irnl_rdt_correction"], capture_output=True, text=True)
17+
18+
# Check if the command executed with sysexit (exit code 2)
19+
assert result.returncode == 2
20+
21+
# Check for expected output
22+
expected_output = "error: the following arguments are required: --beams, --twiss"
23+
assert expected_output in result.stderr
24+
25+
26+
def test_package_execution_fix_v1_1_2():
27+
""" Tests if the package can be run as a module. This failed in <1.1.2 with an import error. """
28+
# Run the package as a module using subprocess
29+
result = subprocess.run([sys.executable, "-m", "irnl_rdt_correction"], capture_output=True, text=True)
30+
31+
# Check if the command executed not with an import error (exit code 1)
32+
assert result.returncode != 1
33+
34+
# Check for a module not found error
35+
not_expected_output = "ModuleNotFoundError: No module named 'utilities'"
36+
assert not_expected_output not in result.stderr
37+
38+
39+
@pytest.mark.parametrize('accel', ('lhc', 'hllhc'))
40+
def test_basic_correction(tmp_path: Path, accel: str):
41+
"""Tests the basic correction functionality and performs some sanity checks.
42+
Same as in tets_standard_correction, but less testcases and called from commandline.
43+
44+
Operates on a pseudo-model so that the corrector values are easily known.
45+
Sanity Checks:
46+
- all correctors found
47+
- correctors have the correct value (as set by errors or zero)
48+
- all corrector circuits present in madx-script
49+
"""
50+
# Parameters -----------------------------------------------------------
51+
order = 4
52+
orientation = ""
53+
54+
correct_ips = (1, 3)
55+
error_value = 2
56+
n_magnets = 4
57+
n_ips = 4
58+
59+
n_correct_ips = len(correct_ips)
60+
n_sides = len("LR")
61+
n_orientation = len(["S", ""])
62+
63+
# Setup ----------------------------------------------------------------
64+
twiss = generate_pseudo_model(accel=accel, n_ips=n_ips, n_magnets=n_magnets)
65+
errors = generate_errortable(index=get_some_magnet_names(n_ips=n_ips, n_magnets=n_magnets))
66+
error_component = f"K{order-1}{orientation}L"
67+
errors[error_component] = error_value
68+
69+
if order % 2: # order is odd -> sides have different sign in rdt
70+
left_hand_magnets = errors.index.str.match(r".*L\d$")
71+
errors.loc[left_hand_magnets, error_component] = errors.loc[left_hand_magnets, error_component] / 2 # so they don't fully compensate
72+
73+
twiss_path = tmp_path / "twiss_input.tfs"
74+
errors_path = tmp_path / "errors_input.tfs"
75+
result_path = tmp_path / "correct"
76+
77+
tfs.write(twiss_path, twiss, save_index="NAME")
78+
tfs.write(errors_path, errors, save_index="NAME")
79+
80+
# Correction -----------------------------------------------------------
81+
subprocess.run([sys.executable, "-m",
82+
"irnl_rdt_correction",
83+
"--accel", accel,
84+
"--twiss", twiss_path,
85+
"--errors", errors_path,
86+
"--beams", "1",
87+
"--output", result_path,
88+
"--feeddown", "0",
89+
"--ips", *[str(ip) for ip in correct_ips],
90+
"--ignore_missing_columns",
91+
"--iterations", "1",
92+
],
93+
capture_output=True, text=True
94+
)
95+
96+
madx_corrections = result_path.with_suffix(".madx").read_text()
97+
df_corrections = tfs.read(result_path.with_suffix(".tfs"))
98+
99+
# Testing --------------------------------------------------------------
100+
# Same as in test_standard_correction - TODO: make function?
101+
# Check output data ---
102+
assert len(list(tmp_path.glob("correct.*"))) == 2
103+
104+
# Check all found correctors ---
105+
if accel == 'lhc':
106+
assert len(df_corrections.index) == (
107+
n_orientation * n_sides * n_correct_ips * len("SO") +
108+
n_sides * n_correct_ips * len("T")
109+
)
110+
111+
if accel == 'hllhc':
112+
assert len(df_corrections.index) == n_orientation * n_sides * n_correct_ips * len("SODT")
113+
114+
# All circuits in madx script ---
115+
for circuit in df_corrections[CIRCUIT]:
116+
assert circuit in madx_corrections
117+
118+
# Check corrector values ---
119+
for test_order in range(3, MAX_N+1):
120+
for test_orientation in ("S", ""):
121+
for ip in correct_ips:
122+
mask = (
123+
(df_corrections[STRENGTH] == f"K{test_order-1}{test_orientation}L") &
124+
(df_corrections[IP] == ip)
125+
)
126+
if (test_order == order) and (test_orientation == orientation):
127+
if order % 2:
128+
corrector_strengths = sum(df_corrections.loc[mask, VALUE])
129+
assert abs(corrector_strengths) < EPS # correctors should be equally distributed
130+
131+
corrector_strengths = -sum(df_corrections.loc[mask, VALUE].abs())
132+
# as beta cancels out (and is 1 anyway)
133+
error_strengths = n_magnets * error_value / 2 # account for partial compensation (from above)
134+
else:
135+
corrector_strengths = sum(df_corrections.loc[mask, VALUE])
136+
assert all(abs(df_corrections.loc[mask, VALUE] - corrector_strengths / n_sides) < EPS)
137+
# as beta cancels out (and is 1 anyway)
138+
error_strengths = (n_sides * n_magnets * error_value)
139+
assert abs(corrector_strengths + error_strengths) < EPS # compensation of RDT
140+
else:
141+
assert all(df_corrections.loc[mask, VALUE] == 0.)

tests/test_dual_optics_correction.py

Lines changed: 2 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -85,7 +85,8 @@ def test_dual_optics(tmp_path: Path):
8585
b1 = beta**(exp_x+exp_y)
8686
b2 = beta2**(exp_x+exp_y)
8787
dual_correction = np.linalg.lstsq(np.array([[b1, b1], [b2, b2]]),
88-
np.array([-b1*error_strengths1, -b2*error_strengths2]))[0]
88+
np.array([-b1*error_strengths1, -b2*error_strengths2]),
89+
rcond=None)[0]
8990

9091
assert all(np.abs(dual_correction) > 0) # just for safety, that there is a solution
9192

tests/test_standard_correction.py

Lines changed: 7 additions & 7 deletions
Original file line numberDiff line numberDiff line change
@@ -228,8 +228,8 @@ class BeamResult:
228228
# here: 2 == sextupole
229229
kl_columns = [f"K{order}{orientation}L" for order in orders for orientation in orientations]
230230

231-
twiss_0.loc[:, kl_columns] = 0
232-
errors_0.loc[:, kl_columns] = 0
231+
twiss_0.loc[:, kl_columns] = 0.0
232+
errors_0.loc[:, kl_columns] = 0.0
233233

234234
# orbit at corrector should not matter, as we don't use feed-down to correct here)
235235
# twiss_0.loc[:, [X, Y]] = 0.5
@@ -243,23 +243,23 @@ class BeamResult:
243243
if magnet_order in ("D", "DS") and is_lhc:
244244
continue
245245
corrector_magnets = twiss_0.index.str.match(rf"MC{magnet_order}X")
246-
twiss_0.loc[corrector_magnets, kl] = 1
246+
twiss_0.loc[corrector_magnets, kl] = 1.0
247247

248248
# Power other magnets and assign errors
249249
non_corrector_magnets = twiss_0.index.str.match(r"M.\.") # M[A-Z]. is created above
250-
twiss_0.loc[non_corrector_magnets, kl_columns] = 1
250+
twiss_0.loc[non_corrector_magnets, kl_columns] = 1.0
251251
twiss_0.loc[non_corrector_magnets, [X, Y]] = 0.5
252252

253253
non_corrector_magnets = errors_0.index.str.match(r"M.\.") # M[A-Z]. is created above
254-
errors_0.loc[non_corrector_magnets, kl_columns] = 1
254+
errors_0.loc[non_corrector_magnets, kl_columns] = 1.0
255255
errors_0.loc[non_corrector_magnets, [f"{DELTA}{X}", f"{DELTA}{Y}"]] = 0.5
256256

257257
# modify the left side, so they don't fully compensate for even (madx)-orders
258258
left_hand_magnets = twiss_0.index.str.match(r".*L\d$")
259-
twiss_0.loc[left_hand_magnets, f"{BETA}{Y}"] = 2 * twiss_0.loc[left_hand_magnets, f"{BETA}{Y}"]
259+
twiss_0.loc[left_hand_magnets, f"{BETA}{Y}"] = 2.0 * twiss_0.loc[left_hand_magnets, f"{BETA}{Y}"]
260260

261261
left_hand_magnets = errors_0.index.str.match(r".*L\d$")
262-
errors_0.loc[left_hand_magnets, kl_columns] = errors_0.loc[left_hand_magnets, kl_columns] / 2
262+
errors_0.loc[left_hand_magnets, kl_columns] = errors_0.loc[left_hand_magnets, kl_columns] / 2.0
263263

264264
# Pre-calculate the integral based on this setup.
265265
integral_left = 1.5 * n_magnets

0 commit comments

Comments
 (0)