Belle II Software  release-08-01-10
test_cluster_position_estimator.py
1 #!/usr/bin/env python3
2 
3 
10 
11 
12 # This steering file steers fills cluster position residuals and pulls with truehits to
13 # test the cluster position estimator payloads from CAF.
14 #
15 # Execute as: basf2 test_cluster_position_estimator.py
16 
17 import math
18 import basf2 as b2
19 import ROOT
20 from ROOT import Belle2
21 
22 # set some random seed
23 b2.set_random_seed(10346)
24 
25 
26 class PXDPositionEstimation(b2.Module):
27  """
28  Histogram the difference position residuals and pulls between clusters and truehits.
29  """
30 
31  def initialize(self):
32  """
33  Create histograms for pulls and residuals
34  """
35 
36 
37  self.nclustersnclusters = 0
38 
39  self.nfound_shapesnfound_shapes = 0
40 
41  self.nfound_offsetnfound_offset = 0
42 
43  # Let's create a root file to store all profiles
44 
45  self.rfilerfile = ROOT.TFile("PXDPositionEstimation.root", "RECREATE")
46  self.rfilerfile.cd()
47 
48 
49  self.hist_map_momentumhist_map_momentum = {}
50 
51  self.hist_map_theta_uhist_map_theta_u = {}
52 
53  self.hist_map_theta_vhist_map_theta_v = {}
54 
55  self.hist_map_clusterchargehist_map_clustercharge = {}
56 
57  # Loop over clusterkinds (=4 means 'all' kinds)
58  for kind in range(5):
59  self.hist_map_momentumhist_map_momentum[kind] = ROOT.TH1F(f"hist_momentum_kind_{kind:d}",
60  f'Particle momentum kind={kind:d}', 5000, 0.0, 10.0)
61  self.hist_map_theta_uhist_map_theta_u[kind] = ROOT.TH1F(f"hist_theta_u_kind_{kind:d}",
62  f'Particle angle thetaU kind={kind:d}', 180, -90, +90)
63  self.hist_map_theta_vhist_map_theta_v[kind] = ROOT.TH1F(f"hist_theta_v_kind_{kind:d}",
64  f'Particle angle thetaV kind={kind:d}', 180, -90, +90)
65  self.hist_map_clusterchargehist_map_clustercharge[kind] = ROOT.TH1F(f"hist_clustercharge_kind_{kind:d}",
66  f'Cluster charge kind={kind:d}', 255, 0.0, 255.0)
67 
68 
69  self.hist_map_residual_uhist_map_residual_u = {}
70 
71  self.hist_map_residual_vhist_map_residual_v = {}
72 
73  self.hist_map_residual_v_specialhist_map_residual_v_special = {}
74 
75  self.hist_map_residual_pull_uhist_map_residual_pull_u = {}
76 
77  self.hist_map_residual_pull_vhist_map_residual_pull_v = {}
78 
79  for kind in range(4):
80  for mode in range(3):
81  self.hist_map_residual_uhist_map_residual_u[(kind, mode)] = ROOT.TH1F(f'hist_map_residual_u_kind_{kind:d}_mode_{mode:d}',
82  f'PXD residual U kind={kind:d} mode={mode:d}',
83  400, -0.007, +0.007)
84  self.hist_map_residual_vhist_map_residual_v[(kind, mode)] = ROOT.TH1F(f'hist_map_residual_v_kind_{kind:d}_mode_{mode:d}',
85  f'PXD residual V kind={kind:d} mode={mode:d}',
86  400, -0.007, +0.007)
87  self.hist_map_residual_pull_uhist_map_residual_pull_u[(kind, mode)] = ROOT.TH1F(f'hist_map_residual_pull_u_kind_{kind:d}_mode_{mode:d}',
88  f'PXD residual pull U kind={kind:d} mode={mode:d}',
89  200, -10, +10)
90  self.hist_map_residual_pull_vhist_map_residual_pull_v[(kind, mode)] = ROOT.TH1F(f'hist_map_residual_pull_v_kind_{kind:d}_mode_{mode:d}',
91  f'PXD residual pull V kind={kind:d} mode={mode:d}',
92  200, -10, +10)
93 
94 
95  self.binlimitsbinlimits = {}
96  self.binlimitsbinlimits[0] = (-90, -30)
97  self.binlimitsbinlimits[1] = (-30, +30)
98  self.binlimitsbinlimits[2] = (+30, +90)
99 
100  for bin in range(3):
101  name = f'hist_map_residual_v_kind_{kind:d}_mode_{mode:d}_special_{bin:d}'
102  title = f'PXD residual V kind={kind:d} mode={mode:d} {self.binlimits[bin][0]:.0f}' + \
103  f'<thetaV<{self.binlimits[bin][1]:.0f}'
104  self.hist_map_residual_v_specialhist_map_residual_v_special[(kind, mode, bin)] = ROOT.TH1F(name, title, 400, -0.007, +0.007)
105 
106  def event(self):
107  """Fill the residual and pull histograms"""
108 
109  # Get truehits
110  truehits = Belle2.PyStoreArray("PXDTrueHits")
111 
112  for truehit in truehits:
113  if isinstance(truehit, Belle2.PXDTrueHit):
114  sensor_info = Belle2.VXD.GeoCache.get(truehit.getSensorID())
115  clusters = truehit.getRelationsFrom("PXDClusters")
116 
117  # now check if we find a cluster
118  for j, cls in enumerate(clusters):
119  # we ignore all clusters where less then 100 electrons come from
120  # our truehit
121  if clusters.weight(j) < 100:
122  continue
123 
124  mom = truehit.getMomentum()
125  tu = mom[0] / mom[2]
126  tv = mom[1] / mom[2]
127  thetaU = math.atan(tu) * 180 / math.pi
128  thetaV = math.atan(tv) * 180 / math.pi
129 
130  # Only look at primary particles -> check if the following is needed
131  # for mcp in truehit.getRelationsFrom("MCParticles"):
132  # if not mcp.hasStatus(1):
133  # reject = True
134 
135  # Get instance of position estimator
137  clusterkind = cls.getKind()
138 
139  # Clusterkinds 0,1,2,3 refer to all cases which can currently
140  # be corrected. Cases where a cluster pixel touches a sensor
141  # edge or contains pixel with varying vPitch are excluded here.
142  if clusterkind <= 3 and mom.Mag() > 0.02:
143 
144  self.nclustersnclusters += 1
145 
146  # Fill momentum and angles for clusterkind
147  self.hist_map_momentumhist_map_momentum[clusterkind].Fill(mom.Mag())
148  self.hist_map_theta_uhist_map_theta_u[clusterkind].Fill(thetaU)
149  self.hist_map_theta_vhist_map_theta_v[clusterkind].Fill(thetaV)
150  self.hist_map_clusterchargehist_map_clustercharge[clusterkind].Fill(cls.getCharge())
151 
152  # Fill clusterkind=4 for all PXD sensors
153  self.hist_map_momentumhist_map_momentum[4].Fill(mom.Mag())
154  self.hist_map_theta_uhist_map_theta_u[4].Fill(thetaU)
155  self.hist_map_theta_vhist_map_theta_v[4].Fill(thetaV)
156  self.hist_map_clusterchargehist_map_clustercharge[4].Fill(cls.getCharge())
157 
158  # Fill the histograms (mode=2)
159  mode = 2
160  pull_u = (truehit.getU() - cls.getU()) / cls.getUSigma()
161  pull_v = (truehit.getV() - cls.getV()) / cls.getVSigma()
162 
163  self.hist_map_residual_uhist_map_residual_u[(clusterkind, mode)].Fill(truehit.getU() - cls.getU())
164  self.hist_map_residual_vhist_map_residual_v[(clusterkind, mode)].Fill(truehit.getV() - cls.getV())
165  self.hist_map_residual_pull_uhist_map_residual_pull_u[(clusterkind, mode)].Fill(pull_u)
166  self.hist_map_residual_pull_vhist_map_residual_pull_v[(clusterkind, mode)].Fill(pull_v)
167 
168  if thetaV >= self.binlimitsbinlimits[0][0] and thetaV < self.binlimitsbinlimits[0][1]:
169  self.hist_map_residual_v_specialhist_map_residual_v_special[(clusterkind, mode, 0)].Fill(truehit.getV() - cls.getV())
170  elif thetaV >= self.binlimitsbinlimits[1][0] and thetaV < self.binlimitsbinlimits[1][1]:
171  self.hist_map_residual_v_specialhist_map_residual_v_special[(clusterkind, mode, 1)].Fill(truehit.getV() - cls.getV())
172  else:
173  self.hist_map_residual_v_specialhist_map_residual_v_special[(clusterkind, mode, 2)].Fill(truehit.getV() - cls.getV())
174 
175  shape_likelyhood = PositionEstimator.getShapeLikelyhood(cls, tu, tv)
176  if shape_likelyhood > 0:
177  self.nfound_shapesnfound_shapes += 1
178 
179  offset = PositionEstimator.getClusterOffset(cls, tu, tv)
180  if offset:
181  # Now, we can safely querry the correction
182  self.nfound_offsetnfound_offset += 1
183 
184  # We need to explicitely add a shift to the offsets
185  # This is not needed when working with PXDRecoHits
186  shiftU = sensor_info.getUCellPosition(cls.getUStart())
187  shiftV = sensor_info.getVCellPosition(cls.getVStart())
188 
189  # Fill the histograms (mode=0)
190  mode = 0
191  pull_u = (truehit.getU() - shiftU - offset.getU()) / (math.sqrt(offset.getUSigma2()))
192  pull_v = (truehit.getV() - shiftV - offset.getV()) / (math.sqrt(offset.getVSigma2()))
193 
194  self.hist_map_residual_uhist_map_residual_u[(clusterkind, mode)].Fill(truehit.getU() - shiftU - offset.getU())
195  self.hist_map_residual_vhist_map_residual_v[(clusterkind, mode)].Fill(truehit.getV() - shiftV - offset.getV())
196  self.hist_map_residual_pull_uhist_map_residual_pull_u[(clusterkind, mode)].Fill(pull_u)
197  self.hist_map_residual_pull_vhist_map_residual_pull_v[(clusterkind, mode)].Fill(pull_v)
198 
199  if thetaV >= self.binlimitsbinlimits[0][0] and thetaV < self.binlimitsbinlimits[0][1]:
200  self.hist_map_residual_v_specialhist_map_residual_v_special[(clusterkind, mode, 0)].Fill(
201  truehit.getV() - shiftV - offset.getV())
202  elif thetaV >= self.binlimitsbinlimits[1][0] and thetaV < self.binlimitsbinlimits[1][1]:
203  self.hist_map_residual_v_specialhist_map_residual_v_special[(clusterkind, mode, 1)].Fill(
204  truehit.getV() - shiftV - offset.getV())
205  else:
206  self.hist_map_residual_v_specialhist_map_residual_v_special[(clusterkind, mode, 2)].Fill(
207  truehit.getV() - shiftV - offset.getV())
208 
209  # Fill the histograms (mode=1)
210  mode = 1
211  self.hist_map_residual_uhist_map_residual_u[(clusterkind, mode)].Fill(truehit.getU() - shiftU - offset.getU())
212  self.hist_map_residual_vhist_map_residual_v[(clusterkind, mode)].Fill(truehit.getV() - shiftV - offset.getV())
213  self.hist_map_residual_pull_uhist_map_residual_pull_u[(clusterkind, mode)].Fill(pull_u)
214  self.hist_map_residual_pull_vhist_map_residual_pull_v[(clusterkind, mode)].Fill(pull_v)
215 
216  if thetaV >= self.binlimitsbinlimits[0][0] and thetaV < self.binlimitsbinlimits[0][1]:
217  self.hist_map_residual_v_specialhist_map_residual_v_special[(clusterkind, mode, 0)].Fill(
218  truehit.getV() - shiftV - offset.getV())
219  elif thetaV >= self.binlimitsbinlimits[1][0] and thetaV < self.binlimitsbinlimits[1][1]:
220  self.hist_map_residual_v_specialhist_map_residual_v_special[(clusterkind, mode, 1)].Fill(
221  truehit.getV() - shiftV - offset.getV())
222  else:
223  self.hist_map_residual_v_specialhist_map_residual_v_special[(clusterkind, mode, 2)].Fill(
224  truehit.getV() - shiftV - offset.getV())
225 
226  else:
227 
228  # Fill the histograms (mode=1)
229  mode = 1
230  pull_u = (truehit.getU() - cls.getU()) / cls.getUSigma()
231  pull_v = (truehit.getV() - cls.getV()) / cls.getVSigma()
232 
233  self.hist_map_residual_uhist_map_residual_u[(clusterkind, mode)].Fill(truehit.getU() - cls.getU())
234  self.hist_map_residual_vhist_map_residual_v[(clusterkind, mode)].Fill(truehit.getV() - cls.getV())
235  self.hist_map_residual_pull_uhist_map_residual_pull_u[(clusterkind, mode)].Fill(pull_u)
236  self.hist_map_residual_pull_vhist_map_residual_pull_v[(clusterkind, mode)].Fill(pull_v)
237 
238  if thetaV >= self.binlimitsbinlimits[0][0] and thetaV < self.binlimitsbinlimits[0][1]:
239  self.hist_map_residual_v_specialhist_map_residual_v_special[(clusterkind, mode, 0)].Fill(truehit.getV() - cls.getV())
240  elif thetaV >= self.binlimitsbinlimits[1][0] and thetaV < self.binlimitsbinlimits[1][1]:
241  self.hist_map_residual_v_specialhist_map_residual_v_special[(clusterkind, mode, 1)].Fill(truehit.getV() - cls.getV())
242  else:
243  self.hist_map_residual_v_specialhist_map_residual_v_special[(clusterkind, mode, 2)].Fill(truehit.getV() - cls.getV())
244 
245  def terminate(self):
246  """
247  Format and write all histograms and plot them
248  """
249 
250  for kind in range(5):
251  self.hist_map_momentumhist_map_momentum[kind].SetLineWidth(2)
252  self.hist_map_momentumhist_map_momentum[kind].SetXTitle('momentum / GeV')
253  self.hist_map_momentumhist_map_momentum[kind].SetYTitle('number of particles')
254 
255  self.hist_map_theta_uhist_map_theta_u[kind].SetLineWidth(2)
256  self.hist_map_theta_uhist_map_theta_u[kind].SetXTitle('thetaU / degree')
257  self.hist_map_theta_uhist_map_theta_u[kind].SetYTitle('number of particles')
258 
259  self.hist_map_theta_vhist_map_theta_v[kind].SetLineWidth(2)
260  self.hist_map_theta_vhist_map_theta_v[kind].SetXTitle('thetaV / degree')
261  self.hist_map_theta_vhist_map_theta_v[kind].SetYTitle('number of particles')
262 
263  self.hist_map_clusterchargehist_map_clustercharge[kind].SetLineWidth(2)
264  self.hist_map_clusterchargehist_map_clustercharge[kind].SetXTitle('cluster charge / ADU')
265  self.hist_map_clusterchargehist_map_clustercharge[kind].SetYTitle('number of particles')
266 
267  for kind in range(4):
268  for mode in range(3):
269  self.hist_map_residual_pull_uhist_map_residual_pull_u[(kind, mode)].SetLineWidth(2)
270  self.hist_map_residual_pull_uhist_map_residual_pull_u[(kind, mode)].SetXTitle('pull u')
271  self.hist_map_residual_pull_uhist_map_residual_pull_u[(kind, mode)].SetYTitle('number of particles')
272 
273  self.hist_map_residual_pull_vhist_map_residual_pull_v[(kind, mode)].SetLineWidth(2)
274  self.hist_map_residual_pull_vhist_map_residual_pull_v[(kind, mode)].SetXTitle('pull v')
275  self.hist_map_residual_pull_vhist_map_residual_pull_v[(kind, mode)].SetYTitle('number of particles')
276 
277  self.hist_map_residual_uhist_map_residual_u[(kind, mode)].SetLineWidth(2)
278  self.hist_map_residual_uhist_map_residual_u[(kind, mode)].SetXTitle('residuals u / cm')
279  self.hist_map_residual_uhist_map_residual_u[(kind, mode)].SetYTitle('number of particles')
280 
281  self.hist_map_residual_vhist_map_residual_v[(kind, mode)].SetLineWidth(2)
282  self.hist_map_residual_vhist_map_residual_v[(kind, mode)].SetXTitle('residuals v / cm')
283  self.hist_map_residual_vhist_map_residual_v[(kind, mode)].SetYTitle('number of particles')
284 
285  for bin in range(3):
286  self.hist_map_residual_v_specialhist_map_residual_v_special[(kind, mode, bin)].SetLineWidth(2)
287  self.hist_map_residual_v_specialhist_map_residual_v_special[(kind, mode, bin)].SetXTitle('residuals v / cm')
288  self.hist_map_residual_v_specialhist_map_residual_v_special[(kind, mode, bin)].SetYTitle('number of particles')
289 
290  hcoverage = ROOT.TH1F("hist_coverage", 'Coverage of corrections', 2, 1, 2)
291  hcoverage.SetBinContent(1, 100.0 * float(self.nfound_offsetnfound_offset / self.nclustersnclusters))
292  hcoverage.SetBinContent(2, 100.0 * float(self.nfound_shapesnfound_shapes / self.nclustersnclusters))
293  hcoverage.SetLineWidth(2)
294  hcoverage.SetYTitle('coverage / %')
295  hcoverage.SetTitle('Coverage of cluster shape corrections')
296 
297  print(f"Coverage of cluster shape corrections is {100.0 * float(self.nfound_offset / self.nclusters):.2f}% ")
298  print(f"Coverage of cluster shape likelyhoods is {100.0 * float(self.nfound_shapes / self.nclusters):.2f}% ")
299 
300  self.rfilerfile.Write()
301  self.rfilerfile.Close()
302 
303 
304 if __name__ == "__main__":
305 
306  import argparse
307  import glob
308  import sys
309 
310  parser = argparse.ArgumentParser(description="Test cluster shape corrections on generic BBbar events")
311  parser.add_argument(
312  '--bglocation',
313  dest='bglocation',
314  default='/home/benjamin/prerelease-01-00-00b-validation/samples',
315  type=str,
316  help='Location of bg overlay files')
317  parser.add_argument('--bkgOverlay', dest='bkgOverlay', action="store_true", help='Perform background overlay')
318  args = parser.parse_args()
319 
320  # Find background overlay files
321  bkgfiles = glob.glob(args.bglocation + '/*.root')
322  if len(bkgfiles) == 0:
323  print('No BG overlay files found')
324  sys.exit()
325 
326  # Now let's create a path to simulate our events.
327  main = b2.create_path()
328  main.add_module("EventInfoSetter", evtNumList=[10000])
329  main.add_module("Gearbox")
330  # We only need the pxd for this
331  main.add_module("Geometry", components=['MagneticField', 'BeamPipe', 'PXD'], useDB=False)
332 
333  # Generate BBbar events
334  main.add_module("EvtGenInput")
335 
336  # Background overlay input
337  if bkgfiles:
338  if args.bkgOverlay:
339  bkginput = b2.register_module('BGOverlayInput')
340  bkginput.param('inputFileNames', bkgfiles)
341  main.add_module(bkginput)
342 
343  main.add_module("FullSim")
344  main.add_module("PXDDigitizer")
345 
346  # Background overlay executor - after digitizer
347  if bkgfiles:
348  if args.bkgOverlay:
349  main.add_module('BGOverlayExecutor', PXDDigitsName='')
350  main.add_module("PXDDigitSorter", digits='')
351 
352  main.add_module("ActivatePXDClusterPositionEstimator")
353  main.add_module("PXDClusterizer")
354 
355  positionestimation = PXDPositionEstimation()
356  main.add_module(positionestimation)
357  main.add_module("Progress")
358 
359  b2.process(main)
360  print(b2.statistics)
Class PXDTrueHit - Records of tracks that either enter or leave the sensitive volume.
Definition: PXDTrueHit.h:31
static PXDClusterPositionEstimator & getInstance()
Main (and only) way to access the PXDClusterPositionEstimator.
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72
static const SensorInfoBase & get(Belle2::VxdID id)
Return a reference to the SensorInfo of a given SensorID.
Definition: GeoCache.h:139
hist_map_clustercharge
Histograms for cluster charge related to particle.
nfound_offset
Counter for clusters where position correction was found in payload.
hist_map_residual_v_special
Histograms for v residuals for smaller thetaV range.
nfound_shapes
Counter for cluster where shape likelyhood was found in payload.