Belle II Software development
BasicSimulationAndReconstruction.py
1#!/usr/bin/env python3
2
3
10
11"""
12<header>
13 <contact>christian.wessel@belle2.org</contact>
14 <input>EvtGenSimNoBkg.root</input>
15 <output>BasicSimulationAndReconstruction.root</output>
16 <description>
17 This validation script checks the basic of our detector simulation (digits)
18 and reconstruction (clusters etc) in a simulation without background
19 </description>
20</header>
21"""
22
23from ROOT import Belle2
24from ROOT import TH1F, TFile, TNamed
25import basf2
26from svd import add_svd_reconstruction
27from pxd import add_pxd_reconstruction
28
29NAME = 'Basic detector simulation and reconstruction'
30CONTACT = 'christian.wessel@belle2.org'
31INPUT_FILE = '../EvtGenSimNoBkg.root'
32OUTPUT_FILE = 'BasicSimulationAndReconstruction.root'
33
34# If False, the script will be called but no events will be processed
35ACTIVE = True
36
37
38def run():
39 """
40 Check the number of hits that are produced by add_simulation and the basic reconstruction in SVD and CDC.
41 """
42 class BasicSimulationAndReconstructionValidation(basf2.Module):
43 """
44 Module to collect information about the number of
45 * PXDDigits
46 * PXDClusters
47 * PXDSpacePoints
48 * SVDShaperDigits
49 * SVDClusters
50 * SVDSpacePoints
51 * CDCHits
52 * CDCHits above threshold
53 * TOPDigits
54 * ARICHDigits
55 * ECLDigits
56 * KLMDigits
57 """
58
59 def __init__(
60 self,
61 name='',
62 contact='',
63 output_file_name=None,
64 ):
65 """Constructor"""
66
67 super().__init__()
68
69 self.validation_name = NAME
70
71 self.contact = CONTACT
72
73 self.output_file_name = OUTPUT_FILE
74
75 def initialize(self):
76 ''' Initialize the Module: book histograms and set descriptions and checks'''
77
78 self.hists = []
79 self.hists.append(TH1F('PXDDigits', 'PXDDigits (no data reduction)',
80 100, 0, 100))
81 self.hists.append(TH1F('PXDClusters', 'PXDClusters (no data reduction)',
82 50, 0, 50))
83 self.hists.append(TH1F('PXDSpacePoints', 'PXDSpacePoints (no data reduction)',
84 50, 0, 50))
85 self.hists.append(TH1F('SVDShaperDigits', 'SVDShaperDigits', 200, 0, 800))
86 self.hists.append(TH1F('SVDClusters', 'SVDClusters', 200, 0, 400))
87 self.hists.append(TH1F('SVDSpacePoints', 'SVDSpacePoints', 200, 0, 600))
88 self.hists.append(TH1F('CDCHits', 'CDCHits (all)', 200, 0, 2000))
89 self.hists.append(TH1F('CDCHitsAboveThreshold', 'CDCHits (above threshold)', 200, 0, 2000))
90 self.hists.append(TH1F('TOPDigits', 'TOPDigits', 200, 0, 800))
91 self.hists.append(TH1F('ARICHDigits', 'ARICHDigits', 200, 0, 200))
92 self.hists.append(TH1F('ECLDigits', 'ECLDigits, m_Amp > 500 (roughly 25 MeV)',
93 100, 0, 100))
94 self.hists.append(TH1F('KLMDigits', 'KLMDigits', 200, 0, 200))
95
96 for h in self.hists:
97 h.SetXTitle(f'Number of {h.GetName()} in event')
98 option = 'shifter'
99 descr = TNamed('Description', 'Number of ' + h.GetName() +
100 ' per event (no BG overlay, just result of add_simulation and basic reconstruction)')
101
102 h.SetYTitle('entries per bin')
103 h.GetListOfFunctions().Add(descr)
104 check = TNamed('Check', 'Distribution must agree with its reference')
105 h.GetListOfFunctions().Add(check)
106 contact = TNamed('Contact', self.contact)
107 h.GetListOfFunctions().Add(contact)
108 options = TNamed('MetaOptions', option)
109 h.GetListOfFunctions().Add(options)
110
111
112 self.nSVDL3 = 768 * 7 * 2 * 2
113
114 def event(self):
115 ''' Event processor: fill histograms '''
116
117 for h in self.hists:
118 digits = Belle2.PyStoreArray(h.GetName())
119 if h.GetName() == 'ECLDigits':
120 n = 0
121 for digit in digits:
122 if digit.getAmp() > 500: # roughly 25 MeV
123 n += 1
124 h.Fill(n)
125 elif h.GetName() == 'CDCHitsAboveThreshold':
126 pass
127 else:
128 h.Fill(digits.getEntries())
129
130 # CDC
131 cdcHits = Belle2.PyStoreArray('CDCHits')
132 nCDC_above_threshold = 0
133 for cdcHit in cdcHits: # count wires of inner/outer layers
134 if cdcHit.getISuperLayer() == 0 and cdcHit.getADCCount() > 15:
135 nCDC_above_threshold += 1
136 if cdcHit.getISuperLayer() != 0 and cdcHit.getADCCount() > 18:
137 nCDC_above_threshold += 1
138 for h in self.hists:
139 if h.GetName() == 'CDCHitsAboveThreshold':
140 h.Fill(nCDC_above_threshold)
141
142 def terminate(self):
143 """ Write histograms to file."""
144
145 tfile = TFile(self.output_file_name, 'recreate')
146 for h in self.hists:
147 h.Write()
148 tfile.Close()
149
150 """ Steering """
151 basf2.set_random_seed(1509)
152
153 path = basf2.create_path()
154
155 path.add_module('RootInput', inputFileName=INPUT_FILE)
156 path.add_module('Gearbox')
157 path.add_module('Geometry')
158 add_svd_reconstruction(path, isROIsimulation=False)
159 add_pxd_reconstruction(path)
160
161 path.add_module(BasicSimulationAndReconstructionValidation())
162
163 path.add_module('Progress')
164
165 print(path)
166 basf2.process(path)
167 print(basf2.statistics)
168
169
170if __name__ == '__main__':
171 if ACTIVE:
172 run()
173 else:
174 print("This validation deactivated and thus basf2 is not executed.\n"
175 "If you want to run this validation, please set the 'ACTIVE' flag above to 'True'.\n"
176 "Exiting.")
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72