Belle II Software  release-08-00-04
trackingInputValidationBkg.py
1 #!/usr/bin/env python3
2 
3 
10 
11 """
12 <header>
13  <contact>software-tracking@belle2.org</contact>
14  <input>EvtGenSim.root</input>
15  <output>TrackingInputValidationBkg.root</output>
16  <description>
17  This module checks the number of input hits for tracking.
18  </description>
19 </header>
20 """
21 
22 from tracking.validation.plot import ValidationPlot
23 from ROOT import Belle2
24 import collections
25 import ROOT
26 import basf2
27 from svd import add_svd_reconstruction
28 from pxd import add_pxd_reconstruction
29 NAME = 'Tracking Input Validation' # not used?
30 CONTACT = 'software-tracking@belle2.org'
31 INPUT_FILE = '../EvtGenSim.root'
32 OUTPUT_FILE = 'TrackingInputValidationBkg.root'
33 N_EVENTS = 1000
34 
35 ACTIVE = True
36 
37 
38 basf2.set_random_seed(1509)
39 
40 
41 class TrackingInputValidation(basf2.Module):
42  """
43  Module to collect information about the number of
44  * PXDDigits
45  * PXDClusters
46  * PXDSpacePoints
47  * SVDShaperDigits
48  * SVDClusters
49  * SVDSpacePoints
50  * CDCHits
51  """
52 
53  def __init__(
54  self,
55  name='',
56  contact='',
57  output_file_name=None,
58  ):
59  """Constructor"""
60 
61  super().__init__()
62 
63  self.validation_namevalidation_name = NAME
64 
65  self.contactcontact = CONTACT
66 
67  self.output_file_nameoutput_file_name = OUTPUT_FILE
68 
69  def initialize(self):
70  """Receive signal at the start of event processing"""
71 
72 
73  self.PXDDigitsPXDDigits = Belle2.PyStoreArray('PXDDigits')
74 
75  self.PXDClustersPXDClusters = Belle2.PyStoreArray('PXDClusters')
76 
77  self.PXDSpacePointsPXDSpacePoints = Belle2.PyStoreArray('PXDSpacePoints')
78 
79 
80  self.SVDShaperDigitsSVDShaperDigits = Belle2.PyStoreArray('SVDShaperDigits')
81 
82  self.SVDClustersSVDClusters = Belle2.PyStoreArray('SVDClusters')
83 
84  self.SVDSpacePointsSVDSpacePoints = Belle2.PyStoreArray('SVDSpacePoints')
85 
86 
87  self.CDCHitsCDCHits = Belle2.PyStoreArray('CDCHits')
88 
89 
90  self.nPXDDigitsnPXDDigits = collections.deque()
91 
92  self.nPXDClustersnPXDClusters = collections.deque()
93 
94  self.nPXDSpacePointsnPXDSpacePoints = collections.deque()
95 
96 
97  self.nSVDShaperDigitsnSVDShaperDigits = collections.deque()
98 
99  self.nSVDClustersnSVDClusters = collections.deque()
100 
101  self.nSVDSpacePointsnSVDSpacePoints = collections.deque()
102 
103 
104  self.nCDCHitsnCDCHits = collections.deque()
105 
106  def event(self):
107  '''Event function'''
108 
109  self.nPXDDigitsnPXDDigits.append(self.PXDDigitsPXDDigits.getEntries())
110  self.nPXDClustersnPXDClusters.append(self.PXDClustersPXDClusters.getEntries())
111  self.nPXDSpacePointsnPXDSpacePoints.append(self.PXDSpacePointsPXDSpacePoints.getEntries())
112  self.nSVDShaperDigitsnSVDShaperDigits.append(self.SVDShaperDigitsSVDShaperDigits.getEntries())
113  self.nSVDClustersnSVDClusters.append(self.SVDClustersSVDClusters.getEntries())
114  self.nSVDSpacePointsnSVDSpacePoints.append(self.SVDSpacePointsSVDSpacePoints.getEntries())
115  self.nCDCHitsnCDCHits.append(self.CDCHitsCDCHits.getEntries())
116 
117  def terminate(self):
118  """Receive signal at the end of event processing"""
119 
120  ''' Saving'''
121  output_tfile = ROOT.TFile(self.output_file_nameoutput_file_name, 'recreate')
122 
123  h_nPXDDigits = ValidationPlot('h_nPXDDigits')
124  h_nPXDDigits.hist(self.nPXDDigitsnPXDDigits, bins=100, lower_bound=0, upper_bound=1000)
125  h_nPXDDigits.contact = self.contactcontact
126  h_nPXDDigits.check = 'Average of 250 +/- 100 is expected'
127  h_nPXDDigits.description = 'Number of PXDDigits after ROI filtering'
128  h_nPXDDigits.title = 'Number of selected PXDDigits per event'
129  h_nPXDDigits.xlabel = 'Number of PXDDigits'
130  h_nPXDDigits.ylabel = ''
131  h_nPXDDigits.write(output_tfile)
132 
133  h_nPXDClusters = ValidationPlot('h_nPXDClusters')
134  h_nPXDClusters.hist(self.nPXDClustersnPXDClusters, bins=100, lower_bound=0, upper_bound=500)
135  h_nPXDClusters.contact = self.contactcontact
136  h_nPXDClusters.check = 'Average of 75 +/- 20 is expected'
137  h_nPXDClusters.description = 'Number of PXDClusters after ROI filtering'
138  h_nPXDClusters.title = 'Number of PXDClusters per event'
139  h_nPXDClusters.xlabel = 'Number of PXDClusters'
140  h_nPXDClusters.ylabel = ''
141  h_nPXDClusters.write(output_tfile)
142 
143  h_nPXDSpacePoints = ValidationPlot('h_nPXDSpacePoints')
144  h_nPXDSpacePoints.hist(self.nPXDSpacePointsnPXDSpacePoints, bins=100, lower_bound=0, upper_bound=500)
145  h_nPXDSpacePoints.contact = self.contactcontact
146  h_nPXDSpacePoints.check = 'Average of 75 +/- 20 is expected'
147  h_nPXDSpacePoints.description = 'Number of PXDSpacePoints after ROI filtering. \
148  Should be the same as the number of PXDClusters'
149  h_nPXDSpacePoints.title = 'Number of PXDSpacePoints per event'
150  h_nPXDSpacePoints.xlabel = 'Number of PXDSpacePoints'
151  h_nPXDSpacePoints.ylabel = ''
152  h_nPXDSpacePoints.write(output_tfile)
153 
154  h_nSVDShaperDigits = ValidationPlot('h_nSVDShaperDigits')
155  h_nSVDShaperDigits.hist(self.nSVDShaperDigitsnSVDShaperDigits, bins=100, lower_bound=0, upper_bound=10000)
156  h_nSVDShaperDigits.contact = self.contactcontact
157  h_nSVDShaperDigits.check = 'First peak at about 3200 +/- 200 is expected, with a second peak round 5500'
158  h_nSVDShaperDigits.description = 'Number of SVDDigits'
159  h_nSVDShaperDigits.title = 'Number of SVDDigits per event'
160  h_nSVDShaperDigits.xlabel = 'Number of SVDDigits'
161  h_nSVDShaperDigits.ylabel = ''
162  h_nSVDShaperDigits.write(output_tfile)
163 
164  h_nSVDClusters = ValidationPlot('h_nSVDClusters')
165  h_nSVDClusters.hist(self.nSVDClustersnSVDClusters, bins=100, lower_bound=0, upper_bound=2000)
166  h_nSVDClusters.contact = self.contactcontact
167  h_nSVDClusters.check = 'Average of 500 +/- 100 is expected'
168  h_nSVDClusters.description = 'Number of SVDClusters'
169  h_nSVDClusters.title = 'Number of SVDClusters per event'
170  h_nSVDClusters.xlabel = 'Number of SVDClusters'
171  h_nSVDClusters.ylabel = ''
172  h_nSVDClusters.write(output_tfile)
173 
174  h_nSVDSpacePoints = ValidationPlot('h_nSVDSpacePoints')
175  h_nSVDSpacePoints.hist(self.nSVDSpacePointsnSVDSpacePoints, bins=100, lower_bound=0, upper_bound=2000)
176  h_nSVDSpacePoints.contact = self.contactcontact
177  h_nSVDSpacePoints.check = 'Average of 200 +/- 50 is expected'
178  h_nSVDSpacePoints.description = 'Number of SVDSpacePoints'
179  h_nSVDSpacePoints.title = 'Number of SVDSpacePoints per event'
180  h_nSVDSpacePoints.xlabel = 'Number of SVDSpacePoints'
181  h_nSVDSpacePoints.ylabel = ''
182  h_nSVDSpacePoints.write(output_tfile)
183 
184  h_nCDCHits = ValidationPlot('h_nCDCHits')
185  h_nCDCHits.hist(self.nCDCHitsnCDCHits, bins=100, lower_bound=0, upper_bound=5000)
186  h_nCDCHits.contact = self.contactcontact
187  h_nCDCHits.check = 'Average of 3000 +/- 200 is expected'
188  h_nCDCHits.description = 'Number of CDCHits'
189  h_nCDCHits.title = 'Number of CDCHits per event'
190  h_nCDCHits.xlabel = 'Number of CDCHits'
191  h_nCDCHits.ylabel = ''
192  h_nCDCHits.write(output_tfile)
193 
194  output_tfile.Close()
195 
196 
197 """ Steering """
198 
199 path = basf2.create_path()
200 
201 path.add_module('RootInput', inputFileName=INPUT_FILE)
202 path.add_module('Gearbox')
203 path.add_module('Geometry')
204 add_svd_reconstruction(path, isROIsimulation=False)
205 add_pxd_reconstruction(path)
206 
207 path.add_module(TrackingInputValidation())
208 
209 path.add_module('Progress')
210 
211 
212 if ACTIVE:
213  print(path)
214  basf2.process(path)
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72
PXDSpacePoints
StoreArray containing the PXDSpacePoints.
def __init__(self, name='', contact='', output_file_name=None)
SVDShaperDigits
StoreArray containing the SVDShaperDigits.
SVDSpacePoints
StoreArray containing the SVDSpacePoints.