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