Belle II Software  release-08-00-04
roiFindingTrackingValidation.py
1 #!/usr/bin/env python3
2 
3 
10 
11 """
12 <header>
13  <contact>software-tracking@belle2.org</contact>
14  <output>ROIFindingValidation.root</output>
15  <description>
16  This module validates the ROI Finding module.
17  </description>
18 </header>
19 """
20 # <input>EvtGenSimNoBkg.root</input>
21 
22 from tracking import add_tracking_for_PXDDataReduction_simulation
23 from tracking.validation.plot import ValidationPlot
24 from ROOT import Belle2
25 import collections
26 import ROOT
27 import basf2
28 from svd import add_svd_reconstruction
29 from simulation import add_simulation
30 NAME = 'ROIFinding' # not used?
31 CONTACT = 'software-tracking@belle2.org'
32 # INPUT_FILE = '../EvtGenSimNoBkg.root' #can't use it because PXDDataReduction in simulated
33 # INPUT_FILE = 'simRootOutput.root' # for debugging purposes
34 OUTPUT_FILE = 'ROIFindingTrackingValidation.root'
35 N_EVENTS = 1000
36 
37 ACTIVE = True
38 
39 
40 basf2.set_random_seed(1509)
41 
42 
44 
45  """Module to collect information about the ROIs (PXDIntercepts & ROIs) and to
46  generate validation plots on the performance of the ROIFinding module."""
47 
48  def __init__(
49  self,
50  name='',
51  contact='',
52  output_file_name=None,
53  nROIs=0
54  ):
55  """Constructor"""
56 
57  super().__init__()
58 
59  self.validation_namevalidation_name = NAME
60 
61  self.contactcontact = CONTACT
62 
63  self.output_file_nameoutput_file_name = OUTPUT_FILE
64 
65 
66  self.nROIsnROIs = 0
67 
68  self.nPXDDigitsnPXDDigits = 0
69 
70  self.nPXDDigitsINnPXDDigitsIN = 0
71 
72  # default binning used for resolution plots over pt
73  # self.resolution_pt_binning = [0.05, 0.1, 0.25, 0.4, 0.6, 1., 1.5, 2., 3., 4.]
74 
75  def initialize(self):
76  """Receive signal at the start of event processing"""
77 
78 
79  self.nROIs_L1nROIs_L1 = collections.deque()
80 
81  self.nROIs_L2nROIs_L2 = collections.deque()
82 
83 
84  self.drfdrf = collections.deque()
85 
86 
87  self.statU_L1statU_L1 = collections.deque()
88 
89  self.statV_L1statV_L1 = collections.deque()
90 
91  self.statU_L2statU_L2 = collections.deque()
92 
93  self.statV_L2statV_L2 = collections.deque()
94 
95 
96  self.theta_phi_L1theta_phi_L1 = collections.deque()
97 
98  def event(self):
99  ''' ROIs quantities'''
100  rois = Belle2.PyStoreArray('ROIs')
101 
102  nL1_ROI = 0
103  nL2_ROI = 0
104  for roi in rois:
105  sensor = roi.getSensorID()
106  if(sensor.getLayerNumber() == 1):
107  nL1_ROI = nL1_ROI + 1
108  else:
109  nL2_ROI = nL2_ROI + 1
110 
111  self.nROIs_L1nROIs_L1.append(nL1_ROI)
112  self.nROIs_L2nROIs_L2.append(nL2_ROI)
113 
114  self.nROIsnROIs = self.nROIsnROIs + rois.getEntries()
115 
116  '''Data Reduction Factor '''
117  tot_pxd = Belle2.PyStoreArray('PXDDigits')
118  in_pxd = Belle2.PyStoreArray('filteredPXDDigits')
119 
120  if len(tot_pxd) > 0:
121  self.drfdrf.append(len(in_pxd) / len(tot_pxd))
122 
123  ''' PXDIntercepts Statistical Error '''
124  inters = Belle2.PyStoreArray('PXDIntercepts')
125  for inter in inters:
126  sensor = Belle2.VxdID(inter.getSensorID())
127  if(sensor.getLayerNumber() == 1):
128  self.statU_L1statU_L1.append(inter.getSigmaU())
129  self.statV_L1statV_L1.append(inter.getSigmaV())
130  else:
131  self.statU_L2statU_L2.append(inter.getSigmaU())
132  self.statV_L2statV_L2.append(inter.getSigmaV())
133 
134  def terminate(self):
135  """Receive signal at the end of event processing"""
136 
137  contact = self.contactcontact
138  basf2.B2RESULT("total nROIs = {}".format(self.nROIsnROIs))
139 
140  ''' Saving'''
141 
142  output_tfile = ROOT.TFile(self.output_file_nameoutput_file_name, 'recreate')
143 
144  h_nROIs_L1 = ValidationPlot('h_nROIs_L1')
145  h_nROIs_L1.hist(self.nROIs_L1nROIs_L1, bins=100, lower_bound=0, upper_bound=100)
146  h_nROIs_L1.contact = contact
147  h_nROIs_L1.check = 'average of 25 +/- 10 is expected'
148  h_nROIs_L1.description = 'ROIs on L1'
149  h_nROIs_L1.title = 'ROIs on L1'
150  h_nROIs_L1.xlabel = 'number of ROIs'
151  h_nROIs_L1.ylabel = ''
152  h_nROIs_L1.write(output_tfile)
153 
154  h_nROIs_L2 = ValidationPlot('h_nROIs_L2')
155  h_nROIs_L2.hist(self.nROIs_L2nROIs_L2, bins=100, lower_bound=0, upper_bound=100)
156  h_nROIs_L2.contact = contact
157  h_nROIs_L2.check = 'average of 25 +/- 10 is expected'
158  h_nROIs_L2.description = 'ROIs on L2'
159  h_nROIs_L2.title = 'ROIs on L2'
160  h_nROIs_L2.xlabel = 'number of ROIs'
161  h_nROIs_L2.ylabel = ''
162  h_nROIs_L2.write(output_tfile)
163 
164  h_drf = ValidationPlot('h_drf')
165  h_drf.hist(self.drfdrf, bins=100, lower_bound=0, upper_bound=1)
166  h_drf.contact = contact
167  h_drf.check = 'with no bkg, the average should be around 70%'
168  h_drf.description = 'Fraction of PXDDDigits inside ROIs'
169  h_drf.title = 'Fraction of PXDDDigits inside ROIs'
170  h_drf.xlabel = 'DRF'
171  h_drf.ylabel = ''
172  h_drf.write(output_tfile)
173 
174  h_statU_L1 = ValidationPlot('h_statU_L1')
175  h_statU_L1.hist(self.statU_L1statU_L1, bins=100, lower_bound=0, upper_bound=0.2)
176  h_statU_L1.contact = contact
177  h_statU_L1.check = 'average should be around 0.021 cm'
178  h_statU_L1.description = 'Statistical Error of the Intercept on L1 planes, along U'
179  h_statU_L1.title = 'L1 Intercept Statistical Uncertainity along U'
180  h_statU_L1.xlabel = 'U stat uncertainity'
181  h_statU_L1.ylabel = ''
182  h_statU_L1.write(output_tfile)
183 
184  h_statV_L1 = ValidationPlot('h_statV_L1')
185  h_statV_L1.hist(self.statV_L1statV_L1, bins=100, lower_bound=0, upper_bound=0.2)
186  h_statV_L1.contact = contact
187  h_statV_L1.check = 'average should be around 0.021 cm'
188  h_statV_L1.description = 'Statistical Error of the Intercept on L1 planes, along V'
189  h_statV_L1.title = 'L1 Intercept Statistical Uncertainity along V'
190  h_statV_L1.xlabel = 'V stat uncertainity (cm)'
191  h_statV_L1.ylabel = ''
192  h_statV_L1.write(output_tfile)
193 
194  h_statU_L2 = ValidationPlot('h_statU_L2')
195  h_statU_L2.hist(self.statU_L2statU_L2, bins=100, lower_bound=0, upper_bound=0.2)
196  h_statU_L2.contact = contact
197  h_statU_L2.check = 'average should be around 0.015 cm'
198  h_statU_L2.description = 'Statistical Error of the Intercept on L2 planes, along U'
199  h_statU_L2.title = 'L2 Intercept Statistical Uncertainity along U'
200  h_statU_L2.xlabel = 'U stat uncertainity (cm)'
201  h_statU_L2.ylabel = ''
202  h_statU_L2.write(output_tfile)
203 
204  h_statV_L2 = ValidationPlot('h_statV_L2')
205  h_statV_L2.hist(self.statV_L2statV_L2, bins=100, lower_bound=0, upper_bound=0.2)
206  h_statV_L2.contact = contact
207  h_statV_L2.check = 'average should be around 0.016 cm'
208  h_statV_L2.description = 'Statistical Error of the Intercept on L2 planes, along V'
209  h_statV_L2.title = 'L2 Intercept Statistical Uncertainity along V'
210  h_statV_L2.xlabel = 'V stat uncertainity (cm)'
211  h_statV_L2.ylabel = ''
212  h_statV_L2.write(output_tfile)
213 
214  output_tfile.Close()
215 
216 
217 """ Steering """
218 
219 path = basf2.create_path()
220 
221 path.add_module('EventInfoSetter', evtNumList=N_EVENTS)
222 path.add_module('EvtGenInput')
223 add_simulation(path, forceSetPXDDataReduction=True, usePXDDataReduction=False)
224 add_svd_reconstruction(path, isROIsimulation=True)
225 
226 # path.add_module('RootInput', inputFileName=INPUT_FILE, entrySequences=["0:{}".format(N_EVENTS-1)])
227 # path.add_module('Gearbox')
228 # path.add_module('Geometry')
229 
230 pxd_unfiltered_digits = 'PXDDigits'
231 pxd_filtered_digits = 'filteredPXDDigits'
232 
233 # SVD tracking
234 svd_reco_tracks = '__ROIsvdRecoTracks'
235 add_tracking_for_PXDDataReduction_simulation(path, ['SVD', 'CDC'], '__ROIsvdClusters') # CDC is not used at the moment!
236 
237 # ROI Finding
238 pxdDataRed = basf2.register_module('PXDROIFinder')
239 param_pxdDataRed = {
240  'recoTrackListName': svd_reco_tracks,
241  'PXDInterceptListName': 'PXDIntercepts',
242  'ROIListName': 'ROIs',
243 }
244 pxdDataRed.param(param_pxdDataRed)
245 path.add_module(pxdDataRed)
246 
247 # Filtering of PXDDigits
248 pxd_digifilter = basf2.register_module('PXDdigiFilter')
249 pxd_digifilter.param('ROIidsName', 'ROIs')
250 pxd_digifilter.param('PXDDigitsName', pxd_unfiltered_digits)
251 pxd_digifilter.param('PXDDigitsInsideROIName', pxd_filtered_digits)
252 pxd_digifilter.param('PXDDigitsOutsideROIName', 'PXDDigitsOutside')
253 path.add_module(pxd_digifilter)
254 
255 ROIValidationPlots = ROIFindingTrackingValidationPlots()
256 path.add_module(ROIValidationPlots)
257 
258 path.add_module('Progress')
259 
260 if ACTIVE:
261  basf2.process(path)
262 
263  print(path)
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
def __init__(self, name='', contact='', output_file_name=None, nROIs=0)
theta_phi_L1
list of the hit angular information in PXD layer 1