Belle II Software development
PXDHitErrors.py
1#!/usr/bin/env python3
2
3
10
11import math
12import basf2 as b2
13
14# Some ROOT tools
15import ROOT
16from ROOT import Belle2
17
18
19class PXDHitErrors(b2.Module):
20
21 """A simple module to check the reconstruction of PXDTrueHits."""
22
23 def __init__(self):
24 """Initialize the module"""
25
26 super().__init__()
27
28 self.file = open('PXDHitErrorOutput.txt', 'w')
29
30 self.vxdid_factors = (8192, 256, 32)
31
32 self.fHisto = ROOT.TFile('PXDPulls.root', 'RECREATE')
33
34 self.h_pull_u = ROOT.TH1F('h_pull_u', 'Pulls in u', 150, -10, 5)
35
36 self.h_pull_v = ROOT.TH1F('h_pull_v', 'Pulls in v', 100, -5, 5)
37
38 def beginRun(self):
39 """ Write legend for file columns """
40
41 self.file.write('vxd.id layer ladder sensor truehit.index cluster.index ')
42 self.file.write('truehit.u truehit.v truehit.time truehit.charge theta.u theta.v ')
43 self.file.write('u v u.error v.error rho charge seed size u.size v.size u.pull v.pull\n')
44
45 def event(self):
46 """Find clusters with a truehit and print some stats."""
47
48 truehits = Belle2.PyStoreArray('PXDTrueHits')
49 # nTruehits = truehits.getEntries()
50 clusters = Belle2.PyStoreArray('PXDClusters')
51 nClusters = clusters.getEntries()
52 # digits = Belle2.PyStoreArray('PXDDigits')
53 # nDigits = digits.getEntries()
54 relClustersToTrueHits = \
55 Belle2.PyRelationArray('PXDClustersToPXDTrueHits')
56 # nClusterRelations = relClustersToTrueHits.getEntries()
57 # relClustersToDigits = Belle2.PyRelationArray('PXDClustersToPXDDigits')
58 # nDigitRelations = relClustersToDigits.getEntries()
59
60 # Start with clusters and use the relation to get the corresponding
61 # digits and truehits.
62 for cluster_index in range(nClusters):
63 cluster = clusters[cluster_index]
64 cluster_truehits = \
65 relClustersToTrueHits.getToIndices(cluster_index)
66
67 # Here we ask only for clusters with exactly one TrueHit.
68 # We don't want combined clusters!
69 if len(cluster_truehits) != 1:
70 continue
71
72 for truehit_index in cluster_truehits:
73 truehit = truehits[truehit_index]
74 # Now let's store some data
75 s = ''
76 # Sesnor identification
77 sensorID = truehit.getRawSensorID()
78 [layer, ladder, sensor] = self.decode(sensorID)
79
80 s_id = f'{sensorID} {layer} {ladder} {sensor} {truehit_index:4d} {cluster_index:4d} '
81 s += s_id
82 # TrueHit information
83 thetaU = math.atan2(truehit.getExitU() - truehit.getEntryU(),
84 0.0075)
85 thetaV = math.atan2(truehit.getExitV() - truehit.getEntryV(),
86 0.0075)
87 s_th = f'{truehit.getU():10.5f} {truehit.getV():10.5f} {truehit.getGlobalTime():10.2f} ' + \
88 f'{truehit.getEnergyDep():10.7f} {thetaU:6.3f} {thetaV:6.3f} '
89 s += s_th
90 # Cluster information
91 cluster_pull_u = 0
92 cluster_pull_v = 0
93 try:
94 cluster_pull_u = (cluster.getU() - truehit.getU()) \
95 / cluster.getUSigma()
96 cluster_pull_v = (cluster.getV() - truehit.getV()) \
97 / cluster.getVSigma()
98 except ZeroDivisionError:
99 if cluster.getUSigma() < 1.0e-8:
100 b2.B2ERROR(f'Zero error in u, clsize {cluster.getUSize()}.')
101 else:
102 b2.B2ERROR(f'Zero error in v, clsize {cluster.getVSize()}.')
103
104 s_cl = \
105 f'{cluster.getU():10.5f} {cluster.getV():10.5f} {cluster.getUSigma():10.5f} {cluster.getVSigma():10.5f} ' + \
106 f'{cluster.getRho():10.4f} {cluster.getCharge():9.1f} {cluster.getSeedCharge():9.1f} ' + \
107 f'{cluster.getSize():5d} {cluster.getUSize():5d} {cluster.getVSize():5d} ' + \
108 f'{cluster_pull_u:10.3f} {cluster_pull_v:10.3f}'
109 s += s_cl
110 # NO DIGITS by now.
111 s += '\n'
112 self.file.write(s)
113 if cluster.getUSize() == 2:
114 self.h_pull_u.Fill(cluster_pull_u)
115 if cluster.getVSize() == 2:
116 self.h_pull_v.Fill(cluster_pull_v)
117
118 def terminate(self):
119 """ Close the output file."""
120
121 self.file.close()
122 # Save histograms to file
123 self.fHisto.cd()
124 self.h_pull_u.Write()
125 self.h_pull_v.Write()
126 self.fHisto.Flush()
127 self.fHisto.Close()
128
129 def decode(self, vxdid):
130 """ Utility to decode sensor IDs """
131
132 result = []
133 for f in self.vxdid_factors:
134 result.append(vxdid / f)
135 vxdid = vxdid % f
136
137 return result
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72
file
Input file object.
Definition: PXDHitErrors.py:28
def decode(self, vxdid)
vxdid_factors
Factors for decoding VXDId's.
Definition: PXDHitErrors.py:30
h_pull_v
Histogram for v pulls.
Definition: PXDHitErrors.py:36
h_pull_u
Histogram for u pulls.
Definition: PXDHitErrors.py:34
fHisto
File to save histograms.
Definition: PXDHitErrors.py:32