Belle II Software development
parameters.py
1
8
9from ROOT import Belle2
10
11
12def cdc_layers(layers=None):
13 if layers is None:
14 layers = [lyr for lyr in range(0, 56)]
15
16 result = []
17
18 wire = 511
19 for layer in layers:
20 for param in [1, 2, 6, 11, 12, 16]:
21 wireid = Belle2.WireID(layer, wire).getEWire()
22 label = Belle2.GlobalLabel()
23 label.construct(Belle2.CDCAlignment.getGlobalUniqueID(), wireid, param)
24 result.append(label.label())
25
26 return result
27
28
29def cdc_wires(layers=None):
30 if layers is None:
31 layers = list(range(0, 56))
32
33 wires_in_layer = [
34 160, 160, 160, 160, 160, 160, 160, 160,
35 160, 160, 160, 160, 160, 160,
36 192, 192, 192, 192, 192, 192,
37 224, 224, 224, 224, 224, 224,
38 256, 256, 256, 256, 256, 256,
39 288, 288, 288, 288, 288, 288,
40 320, 320, 320, 320, 320, 320,
41 352, 352, 352, 352, 352, 352,
42 384, 384, 384, 384, 384, 384]
43
44 result = []
45
46 for layer in layers:
47 for wire in range(0, wires_in_layer[layer]):
48 # Unique id of CDCAlignment db object
50
51 wireid = Belle2.WireID(layer, wire).getEWire()
52 label = Belle2.GlobalLabel()
53
54 label.construct(cdcid, wireid, Belle2.CDCAlignment.wireBwdX)
55 result.append(label.label())
56
57 label.construct(cdcid, wireid, Belle2.CDCAlignment.wireFwdX)
58 result.append(label.label())
59
60 label.construct(cdcid, wireid, Belle2.CDCAlignment.wireBwdY)
61 result.append(label.label())
62
63 label.construct(cdcid, wireid, Belle2.CDCAlignment.wireFwdY)
64 result.append(label.label())
65
66 return result
67
68
69def cdc_t0s():
70 wires_in_layer = [
71 160, 160, 160, 160, 160, 160, 160, 160,
72 160, 160, 160, 160, 160, 160,
73 192, 192, 192, 192, 192, 192,
74 224, 224, 224, 224, 224, 224,
75 256, 256, 256, 256, 256, 256,
76 288, 288, 288, 288, 288, 288,
77 320, 320, 320, 320, 320, 320,
78 352, 352, 352, 352, 352, 352,
79 384, 384, 384, 384, 384, 384]
80
81 result = []
82
83 for layer in range(0, 56):
84 for wire in range(0, wires_in_layer[layer]):
85 label = Belle2.GlobalLabel()
86 label.construct(Belle2.CDCTimeZeros.getGlobalUniqueID(), Belle2.WireID(layer, wire).getEWire(), 0)
87 result.append(label.label())
88
89 return result
90
91
92def vxd_halfshells(pxd=True, svd=True, parameters=None, ying=True, yang=True, pat=True, mat=True):
93 if parameters is None:
94 parameters = [1, 2, 3, 4, 5, 6]
95
96 result = []
97
98 shells = []
99 _ying = Belle2.VxdID(1, 0, 0, 1)
100 _yang = Belle2.VxdID(1, 0, 0, 2)
101 _pat = Belle2.VxdID(3, 0, 0, 1)
102 _mat = Belle2.VxdID(3, 0, 0, 2)
103
104 if pxd:
105 if ying:
106 shells.append(_ying)
107 if yang:
108 shells.append(_yang)
109 if svd:
110 if pat:
111 shells.append(_pat)
112 if mat:
113 shells.append(_mat)
114
115 for vxdid in shells:
116 for param in parameters:
117 vxdid = Belle2.VxdID(vxdid).getID()
118 label = Belle2.GlobalLabel()
119 label.construct(Belle2.VXDAlignment.getGlobalUniqueID(), vxdid, param)
120 result.append(label.label())
121
122 return result
123
124
125def beamspot():
126
127 result = []
128
129 for param in [1, 2, 3]:
130 label = Belle2.GlobalLabel()
131 label.construct(Belle2.BeamSpot.getGlobalUniqueID(), 0, param)
132 result.append(label.label())
133
134 return result
135
136
137def vxd_ladders(layers=None, parameters=None):
138 if layers is None:
139 layers = [1, 2, 3, 4, 5, 6]
140 if parameters is None:
141 parameters = [1, 2, 3, 4, 5, 6]
142
143 result = []
144
145 ladders = [8, 12, 7, 10, 12, 16]
146
147 for layer in layers:
148 for ladder in range(1, ladders[layer - 1] + 1):
149 for ipar in parameters:
150 label = Belle2.GlobalLabel()
151 label.construct(Belle2.VXDAlignment.getGlobalUniqueID(), Belle2.VxdID(layer, ladder, 0).getID(), ipar)
152 result.append(label.label())
153
154 return result
155
156
157def vxd_sensors(layers=None, rigid=True, surface=True, surface2=True, surface3=True, surface4=True, parameters=None):
158 if layers is None:
159 layers = [1, 2, 3, 4, 5, 6]
160
161 params_rigid = [1, 2, 3, 4, 5, 6]
162 params_surface2 = [31, 32, 33]
163 params_surface3 = [41, 42, 43, 44]
164 params_surface4 = [51, 52, 53, 54, 55]
165
166 params = []
167 if rigid:
168 params += params_rigid
169
170 if surface:
171 if surface2:
172 params += params_surface2
173 if surface3:
174 params += params_surface3
175 if surface4:
176 params += params_surface4
177
178 # Allow full override by user
179 if parameters:
180 params = parameters
181
182 result = []
183
184 ladders = [8, 12, 7, 10, 12, 16]
185 sensors = [2, 2, 2, 3, 4, 5]
186
187 for layer in layers:
188 for ladder in range(1, ladders[layer - 1] + 1):
189 for sensor in range(1, sensors[layer - 1] + 1):
190 for param in params:
191 label = Belle2.GlobalLabel()
192 label.construct(Belle2.VXDAlignment.getGlobalUniqueID(), Belle2.VxdID(layer, ladder, sensor).getID(), param)
193 result.append(label.label())
194 return result
195
196
197def vxd():
198 return vxd_sensors() + vxd_ladders() + vxd_halfshells()
199
200
201def pxd():
202 return vxd_sensors(layers=[1, 2]) + vxd_ladders(layers=[1, 2]) + vxd_halfshells(pxd=True, svd=False)
203
204
205def svd():
206 return vxd_sensors(layers=[3, 4, 5, 6]) + vxd_ladders(layers=[3, 4, 5, 6]) + vxd_halfshells(pxd=False, svd=True)
207
208
209def all():
210 # TODO: klm
211 return beamspot() + pxd() + svd() + cdc_layers() + cdc_wires() + cdc_t0s()
212
213
214if __name__ == '__main__':
215 print("Number of available parameters:")
216 print("BeamSpot:", len(beamspot()))
217 print("PXD:", len(pxd()))
218 print("SVD:", len(svd()))
219 print("( VXD:", len(vxd()), ")")
220 assert(len(vxd()) == len(pxd()) + len(svd()))
221 print("CDC layers:", len(cdc_layers()))
222 print("CDC wires:", len(cdc_wires()))
223 print("CDC T0s:", len(cdc_t0s()))
224 print("TOTAL:", len(all()))
static unsigned short getGlobalUniqueID()
Return unique ID of BeamSpot in global Millepede calibration (1)
Definition: BeamSpot.h:91
static unsigned short getGlobalUniqueID()
Get global unique id.
Definition: CDCAlignment.h:109
static unsigned short getGlobalUniqueID()
Get global unique id.
Definition: CDCTimeZeros.h:140
Class to convert to/from global labels for Millepede II to/from detector & parameter identificators.
Definition: GlobalLabel.h:41
static unsigned short getGlobalUniqueID()
Get global unique id.
Definition: VXDAlignment.h:47
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
Class to identify a wire inside the CDC.
Definition: WireID.h:34
Definition: __init__.py:1
Definition: __init__.py:1