ATLAS Offline Software
Loading...
Searching...
No Matches
FPGATrackSimReadTVInputFile Namespace Reference

Functions

 get_regions ()
 print_tree_info (tree, name, check_content=False, max_events=None)
 get_event_number (tree)
 iter_road_hits (road)
 summarize_roads_tracks (header, label, eventLabel, stage, only_nonzero=True)
 check_tree_events (tree, branch_map, max_events=None)
 find_events_with_data (tree, stage)
 main ()

Variables

bool VERBOSE = False

Detailed Description

FPGATrackSimReadTVInputFile.py

Validate that all EDM stages/trees and branches are accessible in the input ROOT file for all regions in the phi slice.
Prints basic info for each TTree and branch.

This script essentially follows the logic used in the bytestreammaker script, but without the bytestreammaker dependencies.

More in https://gitlab.cern.ch/atlas-tdaq-ph2upgrades/atlas-tdaq-eftracking/data-format-tools/bytestreammaker
and https://gitlab.cern.ch/atlas-tdaq-ph2upgrades/atlas-tdaq-eftracking/data-format-tools/eftrackingtestdatautility

Function Documentation

◆ check_tree_events()

FPGATrackSimReadTVInputFile.check_tree_events ( tree,
branch_map,
max_events = None )

Definition at line 128 of file FPGATrackSimReadTVInputFile.py.

128def check_tree_events(tree, branch_map, max_events=None):
129 n_entries = tree.GetEntries()
130 limit = n_entries if max_events is None else min(n_entries, max_events)
131 for i in range(limit):
132 tree.GetEntry(i)
133 evtNum = get_event_number(tree)
134 eventLabel = evtNum if evtNum is not None else f"entry {i}"
135 for bname, info in branch_map.items():
136 if isinstance(info, tuple):
137 label, stage = info
138 else:
139 label, stage = info, None
140 try:
141 header = getattr(tree, bname, None)
142 if header is not None and stage is not None:
143 summarize_roads_tracks(header, label, eventLabel, stage)
144 except Exception as e:
145 raise RuntimeError(f"Event {eventLabel}: {label} [error accessing: {e}]") from e
146
147
#define min(a, b)
Definition cfImp.cxx:40

◆ find_events_with_data()

FPGATrackSimReadTVInputFile.find_events_with_data ( tree,
stage )
Return {event_number: first_hit} for events that have roads+tracks with accessible hits.

Definition at line 148 of file FPGATrackSimReadTVInputFile.py.

148def find_events_with_data(tree, stage):
149 """Return {event_number: first_hit} for events that have roads+tracks with accessible hits."""
150 if not tree:
151 return {}
152 suffix = '2nd' if stage == '2nd' else '1st'
153 events = {}
154 for i in range(tree.GetEntries()):
155 tree.GetEntry(i)
156 evtNum = get_event_number(tree)
157 if evtNum is None:
158 continue
159 header = getattr(tree, 'LogicalEventOutputHeader', None)
160 if header is None:
161 continue
162 roads = getattr(header, f'getFPGATrackSimRoads_{suffix}')()
163 tracks = getattr(header, f'getFPGATrackSimTracks_{suffix}')()
164 if len(roads) == 0 or len(tracks) == 0:
165 continue
166 firstHit = next(
167 (hit for road in roads for hit in iter_road_hits(road)),
168 None
169 )
170 if firstHit is not None:
171 events[evtNum] = firstHit
172 return events
173
174

◆ get_event_number()

FPGATrackSimReadTVInputFile.get_event_number ( tree)
Try to extract an event number from the tree's current entry.

Definition at line 47 of file FPGATrackSimReadTVInputFile.py.

47def get_event_number(tree):
48 """Try to extract an event number from the tree's current entry."""
49 for name in ['LogicalEventSlicedHeader',
50 'LogicalEventStripHeader',
51 'LogicalEventSpacepointHeader',
52 'LogicalEventFirstPixelHeader',
53 'LogicalEventSecondPixelHeader']:
54 header = getattr(tree, name, None)
55 if header is None:
56 continue
57 if hasattr(header, 'event'):
58 try:
59 evt = header.event()
60 if hasattr(evt, 'eventNumber'):
61 return evt.eventNumber()
62 except Exception:
63 pass
64 if hasattr(header, 'getEventNumber'):
65 try:
66 return header.getEventNumber()
67 except Exception:
68 pass
69 return None
70
71

◆ get_regions()

FPGATrackSimReadTVInputFile.get_regions ( )
Return the list of regions in the phi slice.

Definition at line 20 of file FPGATrackSimReadTVInputFile.py.

20def get_regions():
21 """Return the list of regions in the phi slice."""
22 return [34, 98, 162, 226, 290, 354, 418, 482, 546, 610, 674, 738, 802, 866, 930, 994, 1058, 1122, 1186, 1250]
23
24

◆ iter_road_hits()

FPGATrackSimReadTVInputFile.iter_road_hits ( road)
Yield every non-None hit from a single FPGATrackSimRoad.

Definition at line 72 of file FPGATrackSimReadTVInputFile.py.

72def iter_road_hits(road):
73 """Yield every non-None hit from a single FPGATrackSimRoad."""
74 for layerHits in road.getAllHits():
75 for hit in layerHits:
76 if hit is not None:
77 yield hit
78
79
80

◆ main()

FPGATrackSimReadTVInputFile.main ( )

Definition at line 175 of file FPGATrackSimReadTVInputFile.py.

175def main():
176 if len(sys.argv) < 2:
177 raise SystemExit(f"Usage: {sys.argv[0]} <input_root_file>")
178 input_file = sys.argv[1]
179 f = ROOT.TFile.Open(input_file)
180 if not f or f.IsZombie():
181 raise RuntimeError(f"Could not open file: {input_file}")
182 print(f"Opened file: {input_file}\n")
183
184 # Check DataPrep tree
185 tree = f.Get("FPGATrackSimDataPrepTree")
186 print_tree_info(tree, "FPGATrackSimDataPrepTree", check_content=True)
187 if tree:
188 for i in range(tree.GetEntries()):
189 tree.GetEntry(i)
190 for bname in ["LogicalEventInputHeader_PreCluster", "LogicalEventInputHeader_PostCluster"]:
191 try:
192 header = getattr(tree, bname, None)
193 if header is not None:
194 header.event()
195 header.nTowers()
196 except Exception as e:
197 raise RuntimeError(f"Event entry {i}: {bname} [error accessing: {e}]") from e
198
199 test_passed = False
200
201 # Check LogicalEvent trees
202 for region in get_regions():
203 print(f"\n=== Region {region} ===")
204 tree1 = f.Get(f"FPGATrackSimLogicalEventTree_reg{region}")
205 print_tree_info(tree1, f"FPGATrackSimLogicalEventTree_reg{region}", check_content=True)
206 if tree1:
207 check_tree_events(tree1, {
208 "LogicalEventStripHeader": "Strip",
209 "LogicalEventSpacepointHeader": "Spacepoint",
210 "LogicalEventFirstPixelHeader": "FirstPixel",
211 "LogicalEventSecondPixelHeader": "SecondPixel",
212 "LogicalEventOutputHeader": ("Output", "1st")
213 })
214 tree2 = f.Get(f"FPGATrackSimSecondStageTree_reg{region}")
215 print_tree_info(tree2, f"FPGATrackSimSecondStageTree_reg{region}", check_content=True)
216 if tree2:
217 check_tree_events(tree2, {
218 "LogicalEventOutputHeader": ("Output2nd", "2nd"),
219 "LogicalEventSlicedHeader": "Sliced2nd"
220 })
221
222 if not test_passed and tree1 and tree2:
223 events_1st = find_events_with_data(tree1, "1st")
224 events_2nd = find_events_with_data(tree2, "2nd")
225 common_events = sorted(set(events_1st) & set(events_2nd))
226 if common_events:
227 evtNum = common_events[0]
228 hit = events_2nd[evtNum]
229 try:
230 is_real = hit.isReal()
231 except Exception:
232 is_real = False
233 print(
234 f" [PASS] Region {region} event {evtNum}: "
235 f"first hit x={hit.getX()}, y={hit.getY()}, z={hit.getZ()}, isReal={is_real}"
236 )
237 test_passed = True
238
239 f.Close()
240 if not test_passed:
241 raise RuntimeError(
242 "No event found with both 1st and 2nd stage roads/tracks and accessible hits "
243 "in the OutputHeader branch."
244 )
245 print("\nValidation complete.")
246
void print(char *figname, TCanvas *c1)
STL class.
int main()
Definition hello.cxx:18

◆ print_tree_info()

FPGATrackSimReadTVInputFile.print_tree_info ( tree,
name,
check_content = False,
max_events = None )

Definition at line 25 of file FPGATrackSimReadTVInputFile.py.

25def print_tree_info(tree, name, check_content=False, max_events=None):
26 if not tree:
27 print(f" [!] Tree '{name}' not found.")
28 return
29 n_entries = tree.GetEntries()
30 print(f" Tree '{name}': {n_entries} entries")
31 for branch in tree.GetListOfBranches():
32 print(f" Branch: {branch.GetName()} ({branch.GetClassName()})")
33 if check_content and n_entries > 0:
34 limit = n_entries if max_events is None else min(n_entries, max_events)
35 if VERBOSE:
36 print(f" [Checking {limit} events for content...]")
37 for i in range(limit):
38 tree.GetEntry(i)
39 for branch in tree.GetListOfBranches():
40 bname = branch.GetName()
41 try:
42 getattr(tree, bname, None)
43 except Exception as e:
44 raise RuntimeError(f"Event entry {i}: {bname} [error accessing: {e}]") from e
45
46

◆ summarize_roads_tracks()

FPGATrackSimReadTVInputFile.summarize_roads_tracks ( header,
label,
eventLabel,
stage,
only_nonzero = True )
Print summary of roads/tracks for a given stage ('1st' or '2nd').

Definition at line 81 of file FPGATrackSimReadTVInputFile.py.

81def summarize_roads_tracks(header, label, eventLabel, stage, only_nonzero=True):
82 """Print summary of roads/tracks for a given stage ('1st' or '2nd')."""
83 suffix = '2nd' if stage == '2nd' else '1st'
84 roads = getattr(header, f'getFPGATrackSimRoads_{suffix}')()
85 tracks = getattr(header, f'getFPGATrackSimTracks_{suffix}')()
86
87 roadCount = len(roads)
88 roadHitTotal, roadRealTotal = 0, 0
89 for road_idx, road in enumerate(roads):
90 for hit in iter_road_hits(road):
91 roadHitTotal += 1
92 try:
93 if hit.isReal():
94 roadRealTotal += 1
95 except Exception:
96 pass
97 if VERBOSE:
98 try:
99 is_real = hit.isReal()
100 except Exception:
101 is_real = False
102 print(f" Road {road_idx} hit: x={hit.getX():.4f}, y={hit.getY():.4f}, z={hit.getZ():.4f}, isReal={is_real}")
103 if (not only_nonzero) or roadCount > 0:
104 print(f" Event {eventLabel}: {label}: roads={roadCount}, "
105 f"road_hits={roadHitTotal}, road_real_hits={roadRealTotal}")
106
107 trackCount = len(tracks)
108 trackHitTotal, trackRealTotal = 0, 0
109 for track_idx, track in enumerate(tracks):
110 for hit in track.getFPGATrackSimHits():
111 trackHitTotal += 1
112 try:
113 if hit.isReal():
114 trackRealTotal += 1
115 except Exception:
116 pass
117 if VERBOSE:
118 try:
119 is_real = hit.isReal()
120 except Exception:
121 is_real = False
122 print(f" Track {track_idx} hit: x={hit.getX():.4f}, y={hit.getY():.4f}, z={hit.getZ():.4f}, isReal={is_real}")
123 if (not only_nonzero) or trackCount > 0:
124 print(f" Event {eventLabel}: {label}: tracks={trackCount}, "
125 f"track_hits={trackHitTotal}, track_real_hits={trackRealTotal}")
126
127

Variable Documentation

◆ VERBOSE

bool FPGATrackSimReadTVInputFile.VERBOSE = False

Definition at line 17 of file FPGATrackSimReadTVInputFile.py.