FairShip
Loading...
Searching...
No Matches
make_nTuple_SBT.py
Go to the documentation of this file.
1#!/usr/bin/env python
2# SPDX-License-Identifier: LGPL-3.0-or-later
3# SPDX-FileCopyrightText: Copyright CERN for the benefit of the SHiP Collaboration
4
5"""Script to collect muons hitting the SBT (including soft interaction products) to a ROOT file."""
6
7import argparse
8import csv
9import logging
10import os
11
12import ROOT as r
13import shipunit as u
14from tabulate import tabulate
15
16pdg = r.TDatabasePDG.Instance()
17
18
19logging.basicConfig(level=logging.INFO)
20
21
22parser = argparse.ArgumentParser(description=__doc__)
23parser.add_argument("--test", dest="testing_code", help="Run Test", action="store_true")
24parser.add_argument(
25 "-o",
26 "--outputfile",
27 default="muonsProduction_wsoft_SBT.root",
28 help="custom outputfile name",
29)
30parser.add_argument(
31 "-p",
32 "--path",
33 help="path to muon background files",
34 default="/eos/experiment/ship/simulation/bkg/MuonBack_2024helium/8070735",
35)
36args = parser.parse_args()
37
38if args.testing_code:
39 logging.info("test code, output file name overwritten as: muonsProduction_wsoft_SBT_test.root")
40 args.outputfile = "muonsProduction_wsoft_SBT_test.root"
41 selectedmuons = "SelectedMuonsSBT_test.txt"
42else:
43 selectedmuons = "SelectedMuonsSBT.txt"
44
45path = args.path
46logging.info(f"Path to MuonBackground : {path}")
47
48fsel = open(selectedmuons, "w") # noqa: SIM115
49csvwriter = csv.writer(fsel)
50
51output_file = r.TFile.Open(args.outputfile, "recreate")
52output_tree = r.TTree("MuonAndSoftInteractions", "Muon information and soft interaction tracks")
53
54imuondata = r.TVectorD(10) # 10 values: pid, px, py, pz, x, y, z, weight,time_of_hit,nmuons_in_event
55output_tree.Branch("imuondata", imuondata)
56
57track_array = r.TObjArray()
58output_tree.Branch("tracks", track_array)
59
60muon_vetoPoints = r.TClonesArray("vetoPoint")
61output_tree.Branch("muon_vetoPoints", muon_vetoPoints)
62
63muon_UpstreamTaggerPoints = r.TClonesArray("UpstreamTaggerPoint")
64output_tree.Branch("muon_UpstreamTaggerPoints", muon_UpstreamTaggerPoints)
65
66h = {}
67h["PvPt_muon"] = r.TH2F(
68 "PvPt_muon",
69 "The momentum of the muons hitting the SBT(unweighted);P(GeV/c);Pt(GeV/c)",
70 200,
71 0.0,
72 400,
73 40,
74 0.0,
75 20,
76)
77h["n_muon"] = r.TH2I(
78 "n_muon",
79 "Number of muons hitting the SBT per Event;number of mu- per event (unweighted) bin width =1;number of mu+ per event (unweighted) bin width=1",
80 6,
81 0.0,
82 6,
83 6,
84 0.0,
85 6,
86)
87h["n_softtracks"] = r.TH1I("n_softtracks", "Number of soft tracks per muon;;(unweighted)", 200, 0, 2000)
88h["n_sbthits"] = r.TH1I("n_sbthits", "Number of SBT hits per muon;n_sbthits(unweighted);", 900, 0, 900)
89
90
91def printMCTrack(n: int, MCTrack) -> None:
92 """Print MCTrack truth."""
93 mcp = MCTrack[n]
94
95 RED = "\033[91m" # ANSI code Red
96 RESET = "\033[0m" # ANSI code Reset to default
97
98 try:
99 particle_name = pdg.GetParticle(mcp.GetPdgCode()).GetName()
100
101 if particle_name == "mu+" or particle_name == "mu-":
102 particle_name = f"{RED}{particle_name}{RESET} " # Highlight muons in red
103
104 print(
105 " %6s %-10s %10i %6.3F %6.3F %7.3F %7.3F %7.3F %7.3F %6s %10.3F %28s"
106 % (
107 n,
108 particle_name,
109 mcp.GetPdgCode(),
110 mcp.GetPx() / u.GeV,
111 mcp.GetPy() / u.GeV,
112 mcp.GetPz() / u.GeV,
113 mcp.GetStartX() / u.m,
114 mcp.GetStartY() / u.m,
115 mcp.GetStartZ() / u.m,
116 mcp.GetMotherId(),
117 mcp.GetWeight(),
118 mcp.GetProcName().Data(),
119 )
120 )
121 except Exception:
122 print(
123 " %6s %-10s %10i %6.3F %6.3F %7.3F %7.3F %7.3F %7.3F %6s %10.3F %28s"
124 % (
125 n,
126 "----",
127 mcp.GetPdgCode(),
128 mcp.GetPx() / u.GeV,
129 mcp.GetPy() / u.GeV,
130 mcp.GetPz() / u.GeV,
131 mcp.GetStartX() / u.m,
132 mcp.GetStartY() / u.m,
133 mcp.GetStartZ() / u.m,
134 mcp.GetMotherId(),
135 mcp.GetWeight(),
136 mcp.GetProcName().Data(),
137 )
138 )
139
140
141def dump(event, pcut: float = 0, print_whole_event: bool = True) -> None:
142 """Dump the whole event."""
143 if print_whole_event:
144 print(
145 "\n %6s %-10s %10s %6s %6s %7s %7s %7s %7s %6s %10s %18s"
146 % (
147 "#",
148 "particle",
149 "pid",
150 "px",
151 "py",
152 "pz",
153 "vx",
154 "vy",
155 "vz",
156 "mid",
157 "w",
158 "Process",
159 )
160 )
161 print(
162 " %6s %10s %10s %6s %6s %7s %7s %7s %7s %6s %10s %18s\n "
163 % (
164 " ",
165 "--------",
166 "---",
167 "--",
168 "--",
169 "--",
170 "--",
171 "--",
172 "--",
173 "---",
174 "---",
175 "-------",
176 )
177 )
178 n = -1
179 for mcp in event.MCTrack:
180 n += 1
181 if mcp.GetP() / u.GeV < pcut:
182 continue
183 if print_whole_event:
184 printMCTrack(n, event.MCTrack)
185
186 return
187
188
189global_event_nr = 0
190processed_events = {}
191P_threshold = 3
192
193headers = [
194 f"nMuons in event>{P_threshold}GeV",
195 "Muon PID",
196 "Momentum[GeV/c]",
197 "x[m]",
198 "y[m]",
199 "z[m]",
200 "t_muon [ns]",
201 "nSoft Tracks",
202 "nSBT Hits(muon)",
203 "nUBT Hits(muon)",
204 "Weight_muon",
205]
206
207csvwriter.writerow(headers)
208
209for inputFolder in os.listdir(path):
210 if not os.path.isdir(os.path.join(path, inputFolder)):
211 continue
212
213 if args.testing_code and global_event_nr >= 100000:
214 break
215
216 logging.info(f"Processing folder: {inputFolder}")
217
218 f = None
219 try:
220 f = r.TFile.Open(
221 os.path.join(path, inputFolder, "ship.conical.MuonBack-TGeant4.root"),
222 "read",
223 )
224 tree = f.cbmsim
225 except Exception as e:
226 logging.debug(f"Error :{e}")
227
228 if f:
229 f.Close()
230 continue
231
232 for eventNr, event in enumerate(tree):
233 global_event_nr += 1
234
235 muon_table = []
236
237 nmu = {"mu+": 0, "mu-": 0}
238
239 muon_ids = []
240 muon_hits = {}
241 # Collect track IDs of muons which hit the SBT
242 for hit in event.vetoPoint:
243 detID = hit.GetDetectorID()
244 pid = hit.PdgCode()
245 track_id = hit.GetTrackID()
246 P = r.TMath.Sqrt(hit.GetPx() ** 2 + hit.GetPy() ** 2 + hit.GetPz() ** 2)
247
248 if 1000 < detID < 999999 and abs(pid) == 13 and P_threshold / u.GeV < P:
249 particle_name = pdg.GetParticle(hit.PdgCode()).GetName()
250 if track_id not in muon_ids:
251 muon_ids.append(track_id)
252 muon_hits[track_id] = 0
253 nmu[particle_name] += 1
254 muon_hits[track_id] += 1
255
256 if not len(muon_ids):
257 continue
258
259 logging.debug(f"\n\nEVENT ID:{global_event_nr}")
260 logging.debug(f"Muon Track Available:{muon_hits.keys()}\n Number of SBT hits within:{muon_hits}")
261
262 for muon_ in muon_ids:
263 imuondata.Zero()
264
265 # saving soft tracks
266 track_array.Clear()
267
268 for track in event.MCTrack:
269 if track.GetMotherId() == muon_ and track.GetProcName().Data() != "Muon nuclear interaction":
270 track_array.Add(track)
271
272 # saving the muon info
273
274 muon_UpstreamTaggerPoints.Clear()
275
276 ubt_index = 0
277
278 for hit in event.UpstreamTaggerPoint:
279 track_id = hit.GetTrackID()
280
281 if track_id != muon_:
282 continue
283
284 if muon_UpstreamTaggerPoints.GetSize() == ubt_index:
285 muon_UpstreamTaggerPoints.Expand(ubt_index + 1)
286 muon_UpstreamTaggerPoints[ubt_index] = hit
287
288 ubt_index += 1
289
290 index = 0
291
292 muon_vetoPoints.Clear()
293
294 for hit in event.vetoPoint:
295 detID = hit.GetDetectorID()
296 pid = hit.PdgCode()
297 track_id = hit.GetTrackID()
298
299 P = r.TMath.Sqrt(hit.GetPx() ** 2 + hit.GetPy() ** 2 + hit.GetPz() ** 2)
300 Pt = r.TMath.Sqrt(hit.GetPx() ** 2 + hit.GetPy() ** 2)
301
302 if 1000 < detID < 999999 and track_id == muon_ and P_threshold / u.GeV < P:
303 if global_event_nr not in processed_events:
304 processed_events[global_event_nr] = []
305
306 if muon_vetoPoints.GetSize() == index:
307 muon_vetoPoints.Expand(index + 1)
308 muon_vetoPoints[index] = hit
309
310 index += 1
311
312 weight = event.MCTrack[track_id].GetWeight()
313
314 if track_id not in processed_events[global_event_nr]: # only save the info of first SBT hit
315 processed_events[global_event_nr].append(track_id)
316
317 imuondata[0] = float(pid)
318 imuondata[1] = float(hit.GetPx() / u.GeV)
319 imuondata[2] = float(hit.GetPy() / u.GeV)
320 imuondata[3] = float(hit.GetPz() / u.GeV)
321 imuondata[4] = float(hit.GetX() / u.m)
322 imuondata[5] = float(hit.GetY() / u.m)
323 imuondata[6] = float(hit.GetZ() / u.m)
324 imuondata[7] = float(weight)
325 imuondata[8] = float(hit.GetTime())
326 imuondata[9] = len(muon_ids)
327 muon_table.append(
328 [
329 track_id,
330 len(muon_ids),
331 imuondata[0],
332 P,
333 imuondata[4],
334 imuondata[5],
335 imuondata[6],
336 imuondata[8],
337 len(track_array),
338 len(muon_vetoPoints),
339 len(muon_UpstreamTaggerPoints),
340 imuondata[7],
341 ]
342 )
343
344 h["PvPt_muon"].Fill(P, Pt)
345 h["n_softtracks"].Fill(len(track_array))
346
347 muon_table[-1][-2] = len(muon_vetoPoints)
348
349 h["n_sbthits"].Fill(len(muon_vetoPoints))
350
351 output_tree.Fill()
352
353 h["n_muon"].Fill(nmu["mu-"], nmu["mu+"])
354 csvwriter.writerows(row[1:] for row in muon_table)
355 # dump(event)
356 logging.debug(f"Muon Summary:\n{tabulate(muon_table, headers=headers, tablefmt='grid')}\n\n")
357
358 f.Close()
359
360
361total_muons = sum(len(values) for values in processed_events.values())
362print(f"nMuons saved: {total_muons}, File: {args.outputfile}")
363
364output_file.cd()
365output_tree.Write()
366for histname in h:
367 h[histname].Write()
368output_file.Close()
369fsel.close()
370
371print(
372 "------------------------------------------------------file saved, inspecting",
373 args.outputfile,
374 "now----------------------------------------------------------------",
375)
376
377event_data = []
378with r.TFile.Open(args.outputfile, "read") as file:
379 try:
380 tree = file.MuonAndSoftInteractions
381 except Exception as e:
382 print(f"Error: {e}")
383 exit(1)
384
385 print(f"Processing tree: {tree.GetName()}")
386 print(f"Total number of entries: {tree.GetEntries()}")
387
388 for event in tree:
389 imuondata = event.imuondata
390 pid = imuondata[0]
391 px = imuondata[1]
392 py = imuondata[2]
393 pz = imuondata[3]
394 x = imuondata[4]
395 y = imuondata[5]
396 z = imuondata[6]
397
398 weight = imuondata[7]
399 time_hit = imuondata[8]
400 nmuons = imuondata[9]
401
402 num_tracks = len(event.tracks)
403 num_muonhits = len(event.muon_vetoPoints)
404 num_ubthits = len(event.muon_UpstreamTaggerPoints)
405
406 P = r.TMath.Sqrt(px**2 + py**2 + pz**2)
407
408 event_data.append(
409 [
410 nmuons,
411 pid,
412 P,
413 x,
414 y,
415 z,
416 time_hit,
417 num_tracks,
418 num_muonhits,
419 num_ubthits,
420 weight,
421 ]
422 )
423
424print(tabulate(event_data, headers=headers, tablefmt="grid"))
None printMCTrack(int n, MCTrack)
None dump(event, float pcut=0, bool print_whole_event=True)
int open(const char *, int)
Opens a file descriptor.