Belle II Software development
vxdalignMon.py
1#!/usr/bin/env python3
2
3
10
11"""
12Generate alignment Monitoring Objects for Mirabelle
13"""
14
15
16import os
17import re
18import ROOT
19from ROOT import Belle2
20import numpy as np
21from glob import glob
22import sys
23
24# extract alignment parameter values from the VXDAlignment payload
25
26
27def get_means(payload):
28 # sensor IDs for PXD and SVD
29 pxd_shells = [Belle2.VxdID(1, 0, 0, 1), Belle2.VxdID(1, 0, 0, 2)]
30 svd_shells = [Belle2.VxdID(3, 0, 0, 1), Belle2.VxdID(3, 0, 0, 2)]
31 means = {}
32 for par in range(1, 7):
33 pxd_values = []
34 svd_values = []
35 vxd_values = []
36 for shell in pxd_shells:
37 v = payload.get(shell.getID(), par)
38 if v != 0.:
39 v = v * 1e3 if par in [4, 5, 6] else v * 1e4
40 pxd_values.append(v)
41 vxd_values.append(v)
42 for shell in svd_shells:
43 v = payload.get(shell.getID(), par)
44 if v != 0.:
45 v = v * 1e3 if par in [4, 5, 6] else v * 1e4
46 svd_values.append(v)
47 vxd_values.append(v)
48 param_name = {1: "x", 2: "y", 3: "z", 4: "alpha", 5: "beta", 6: "gamma"}[par]
49 means[f"pxd_{param_name}"] = np.mean(pxd_values) if pxd_values else 0.0
50 means[f"svd_{param_name}"] = np.mean(svd_values) if svd_values else 0.0
51 means[f"vxd_{param_name}"] = np.mean(vxd_values) if vxd_values else 0.0
52 return means
53
54
55def create_output_file(exp, run, means, out_dir):
56 # create a root file with monitoringObject and DQMFileMetaData for a specific run
57 mon = ROOT.Belle2.MonitoringObject("vxd_alignment")
58 for name, value in means.items():
59 mon.setVariable(name, float(value))
60 meta = ROOT.Belle2.DQMFileMetaData()
61 meta.setExperimentRun(exp, run)
62 if not os.path.exists(out_dir):
63 os.makedirs(out_dir)
64 out_file_path = os.path.join(out_dir, f"mon_exp{exp:03d}_run{run:06d}.root")
65 out_file = ROOT.TFile(out_file_path, "RECREATE")
66 meta.Write()
67 mon.Write()
68 out_file.Close()
69 print(f"Created {out_file_path}")
70
71
72if len(sys.argv) != 3:
73 print("Usage: python merged_code.py <input_path> <output_dir>")
74 sys.exit(1)
75
76input_path = sys.argv[1] # Directory with database.txt and payloads from DDB
77output_dir = sys.argv[2]
78
79# Read database.txt to map revision IDs to experiment/run ranges
80database_file = os.path.join(input_path, "database.txt")
81if not os.path.exists(database_file):
82 print(f"Error: {database_file} not found.")
83 sys.exit(1)
84
85with open(database_file, "r") as db_file:
86 database_lines = db_file.readlines()
87
88# extract run range information from database.txt
89database_dict = {}
90for line in database_lines:
91 match = re.match(r"dbstore/VXDAlignment\s+([A-Za-z0-9]{4})\s+(\d+),(\d+),(\d+),(\d+)", line.strip())
92 if match:
93 rev = match.group(1)
94 exp1 = int(match.group(2))
95 start_run = int(match.group(3))
96 exp2 = int(match.group(4))
97 end_run = int(match.group(5))
98 database_dict[rev] = (exp1, start_run, exp2, end_run)
99
100# Process each .root file
101root_files = glob(os.path.join(input_path, "dbstore_VXDAlignment_rev_*.root"))
102for root_file in root_files:
103 match = re.search(r"dbstore_VXDAlignment_rev_([A-Za-z0-9]{4})\.root", root_file)
104 if not match:
105 print(f"Couldn't find a rev from file name: {root_file}")
106 continue
107 rev = match.group(1)
108 if rev not in database_dict:
109 print(f"No database entry found for rev: {rev}")
110 continue
111 # experiment and run range info from database
112 exp1, start_run, exp2, end_run = database_dict[rev]
113 if exp1 != exp2:
114 print(f"Warning: exp1 ({exp1}) differs from exp2 ({exp2}) for rev {rev}. Using exp1.")
115 exp = exp1
116
117 f = ROOT.TFile.Open(root_file)
118 if not f or not f.IsOpen():
119 print(f"error: unable to open ROOT file: {root_file}")
120 continue
121 payload = f.Get("VXDAlignment")
122 if not payload:
123 print(f"error: VXDAlignment not found in {root_file}")
124 f.Close()
125 continue
126
127 # Calculate mean alignment parameters
128 means = get_means(payload)
129 f.Close()
130
131 # Generate output files for each run
132 for run in range(start_run, end_run + 1):
133 create_output_file(exp, run, means, output_dir)
134
135print(f" Done! all ROOT files created in {output_dir}")
Class to uniquely identify a any structure of the PXD and SVD.
Definition VxdID.h:32