Belle II Software development
Ntuple Class Reference
Inheritance diagram for Ntuple:

Public Member Functions

def initialize (self)
 
def sortPeaks (self, unsortedPeaks)
 
def pdfHistogram (self, pdf)
 
def event (self)
 
def terminate (self)
 

Public Attributes

 trk_selector
 track selector - selection of muons from ee->mumu events
 
 file
 file object
 
 bunchOffset
 histogram of bunch offset
 
 h_cth_vs_p
 histogram cos(theta) vs, momentum
 
 h_momentum
 histogram of momentum
 
 h_cth
 histogram of cos(theta)
 
 h_phi
 histogram of local phi
 
 h_z
 histogram of local z
 
 h_x
 histogram of local x
 
 h_charge
 histogram of charge
 
 h_poca_xy
 histogram of POCA x-y
 
 h_poca_z
 histogram of POCA z
 
 h_hitsCDC
 histogram of number of hits in CDC
 
 h_Ecms
 histogram of Ecms
 
 tree
 tree object
 
 data
 data structure
 

Static Public Attributes

int nhisto = 0
 histogram counter
 

Detailed Description

 Makes a flat ntuple 

Definition at line 99 of file cdst_timeResoNtuple.py.

Member Function Documentation

◆ event()

def event (   self)
 event processing 

Definition at line 203 of file cdst_timeResoNtuple.py.

203 def event(self):
204 ''' event processing '''
205
206 recBunch = Belle2.PyStoreObj('TOPRecBunch')
207 if not recBunch:
208 b2.B2ERROR('no TOPRecBunch')
209 return
210 if not recBunch.isReconstructed():
211 return
212 self.bunchOffset.Fill(recBunch.getCurrentOffset())
213
214 evtMetaData = Belle2.PyStoreObj('EventMetaData')
215 self.data.run = evtMetaData.getRun()
216 self.data.offset = recBunch.getCurrentOffset()
217 self.data.usedTrk = recBunch.getUsedTracks()
218
219 for track in Belle2.PyStoreArray('Tracks'):
220 trk = Belle2.TOP.TOPTrack(track)
221 if not trk.isValid():
222 continue
223 if not self.trk_selector.isSelected(trk):
224 continue
225 pdfs = track.getRelated('TOPPDFCollections')
226 if not pdfs:
227 continue
228 self.data.slot = pdfs.getModuleID()
229 momentum = pdfs.getAssociatedLocalMomentum()
230 position = pdfs.getAssociatedLocalHit()
231 self.data.p = momentum.R()
232 self.data.cth = math.cos(momentum.Theta())
233 self.data.phi = momentum.Phi()
234 self.data.z = position.Z()
235 self.data.x = position.X()
236 self.data.tof = trk.getTOF(Belle2.Const.muon)
237 try:
238 tfit = track.getTrackFitResultWithClosestMass(Belle2.Const.muon)
239 except BaseException:
240 b2.B2ERROR("No trackFitResult available")
241 continue
242 self.data.charge = tfit.getChargeSign()
243 pocaPosition = tfit.getPosition()
244 self.data.poca_x = pocaPosition.X()
245 self.data.poca_y = pocaPosition.Y()
246 self.data.poca_z = pocaPosition.Z()
247 self.data.hitsCDC = tfit.getHitPatternCDC().getNHits()
248 self.data.Ecms = self.trk_selector.getCMSEnergy()
249 try:
250 extHit = trk.getExtHit()
251 timeZero = extHit.getRelated('TOPTimeZeros')
252 self.data.rec_t0 = timeZero.getTime()
253 self.data.valid_t0 = timeZero.isValid()
254 except BaseException:
255 self.data.rec_t0 = 0
256 self.data.valid_t0 = 0
257
258 self.h_cth_vs_p.Fill(self.data.p, self.data.cth)
259 self.h_momentum.Fill(self.data.p)
260 self.h_cth.Fill(self.data.cth)
261 self.h_phi.Fill(self.data.phi)
262 self.h_z.Fill(self.data.z)
263 self.h_x.Fill(self.data.x)
264 self.h_charge.Fill(self.data.charge)
265 self.h_poca_xy.Fill(self.data.poca_x, self.data.poca_y)
266 self.h_poca_z.Fill(self.data.poca_z)
267 self.h_hitsCDC.Fill(self.data.hitsCDC)
268 self.h_Ecms.Fill(self.data.Ecms)
269 try:
270 pdf = pdfs.getHypothesisPDF(13)
271 except BaseException:
272 b2.B2ERROR("No PDF available for PDG = 13")
273 continue
274 self.data.itrk += 1
275 self.pdfHistogram(pdf)
276 emi = trk.getEmissionPoint()
277 x0 = emi.position.X()
278 z0 = emi.position.Y()
279 self.data.nfot = 0
280 for digit in Belle2.PyStoreArray('TOPDigits'):
281 if digit.getModuleID() == self.data.slot and digit.getHitQuality() == 1:
282 k = self.data.nfot
283 if k >= 1000:
284 continue
285 self.data.nfot += 1
286 self.data.channel[k] = digit.getChannel()
287 self.data.pixel[k] = digit.getPixelID()
288 self.data.time[k] = digit.getTime()
289 self.data.timeErr[k] = digit.getTimeError()
290 self.data.pulseHeight[k] = digit.getPulseHeight()
291 self.data.pulseWidth[k] = digit.getPulseWidth()
292 self.data.sample[k] = digit.getModulo256Sample()
293 self.data.status[k] = digit.getStatus()
294 peaks = pdf[digit.getPixelID() - 1]
295 if peaks.empty():
296 self.data.t0[k] = 0
297 self.data.wid0[k] = 0
298 self.data.t1[k] = 0
299 else:
300 sorted_peaks = self.sortPeaks(peaks)
301 self.data.t0[k] = sorted_peaks[0].mean
302 self.data.wid0[k] = sorted_peaks[0].width
303 self.data.t1[k] = self.data.t0[k] + 100
304 if peaks.size() > 1:
305 self.data.t1[k] = sorted_peaks[1].mean
306 self.data.t_pdf[k] = 0
307 self.data.wid[k] = 0
308 self.data.fic[k] = 0
309 self.data.type[k] = 0
310 self.data.nx[k] = 0
311 self.data.ny[k] = 0
312 self.data.nxm[k] = 0
313 self.data.nym[k] = 0
314 self.data.nxe[k] = 0
315 self.data.nye[k] = 0
316 self.data.nys[k] = 0
317 self.data.xd[k] = 0
318 self.data.yd[k] = 0
319 self.data.xm[k] = 0
320 self.data.kx[k] = 0
321 self.data.ky[k] = 0
322 self.data.alpha[k] = 0
323 assocPDF = digit.getRelated('TOPAssociatedPDFs')
324 if assocPDF:
325 pik = assocPDF.getSinglePeak()
326 if pik:
327 self.data.t_pdf[k] = pik.position
328 self.data.wid[k] = pik.width
329 self.data.fic[k] = math.degrees(pik.fic)
330 self.data.type[k] = pik.type
331 self.data.nx[k] = pik.nx
332 self.data.ny[k] = pik.ny
333 self.data.nxm[k] = pik.nxm
334 self.data.nym[k] = pik.nym
335 self.data.nxe[k] = pik.nxe
336 self.data.nye[k] = pik.nye
337 if pik.kyd > 0:
338 self.data.nys[k] = int(pik.nye / 2)
339 else:
340 self.data.nys[k] = int((pik.nye + 1) / 2)
341 self.data.xd[k] = pik.xd
342 self.data.yd[k] = pik.yd
343 if pik.kze > 0:
344 self.data.xm[k] = x0 + pik.kxe / pik.kze * (130.0 - z0)
345 self.data.kx[k] = pik.kxe
346 self.data.ky[k] = pik.kye
347 self.data.alpha[k] = math.degrees(math.acos(abs(pik.kzd)))
348
349 self.tree.Fill()
350
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72
a (simplified) python wrapper for StoreObjPtr.
Definition: PyStoreObj.h:67
Reconstructed track at TOP.
Definition: TOPTrack.h:39

◆ initialize()

def initialize (   self)
 initialize: open root file, construct ntuple 

Definition at line 105 of file cdst_timeResoNtuple.py.

105 def initialize(self):
106 ''' initialize: open root file, construct ntuple '''
107
108 evtMetaData = Belle2.PyStoreObj('EventMetaData')
109 expNo = evtMetaData.obj().getExperiment()
110 runNo = evtMetaData.obj().getRun()
111 exp_run = '-e' + f'{expNo:04d}' + '-r' + f'{runNo:05d}'
112 outName = 'timeResoNtuple' + exp_run + '.root'
113
114
115 self.trk_selector = Belle2.TOP.TrackSelector("dimuon")
116 self.trk_selector.setDeltaEcms(0.1)
117 self.trk_selector.setCutOnPOCA(2.0, 4.0)
118 self.trk_selector.setCutOnLocalZ(-130.0, 130.0)
119
120
121 self.file = ROOT.TFile(outName, 'recreate')
122
123
124 self.bunchOffset = TH1F("bunchOffset", "bunch offset", 100, -1.0, 1.0)
125 self.bunchOffset.SetXTitle("bunch offset [ns]")
126
127 self.h_cth_vs_p = TH2F("cth_vs_p", "local cos #theta vs. p", 100, 0.0, 10.0,
128 100, -1.0, 1.0)
129 self.h_cth_vs_p.SetXTitle("p [GeV/c]")
130 self.h_cth_vs_p.SetYTitle("cos #theta")
131
132 self.h_momentum = TH1F("momentum", "momentum", 100, 0.0, 10.0)
133 self.h_momentum.SetXTitle("p [GeV/c]")
134
135 self.h_cth = TH1F("cos_theta", "local cos #theta", 100, -1.0, 1.0)
136 self.h_cth.SetXTitle("cos #theta")
137
138 self.h_phi = TH1F("local_phi", "local phi", 100, -math.pi, math.pi)
139 self.h_phi.SetXTitle("local #phi")
140
141 self.h_z = TH1F("local_z", "local z", 100, -140.0, 140.0)
142 self.h_z.SetXTitle("local z [cm]")
143
144 self.h_x = TH1F("local_x", "local x", 100, -23.0, 23.0)
145 self.h_x.SetXTitle("local x [cm]")
146
147 self.h_charge = TH1F("charge", "charge", 3, -1.5, 1.5)
148 self.h_charge.SetXTitle("charge")
149
150 self.h_poca_xy = TH2F("poca_xy", "POCA distribution in x-y", 100, -1.0, 1.0,
151 100, -1.0, 1.0)
152 self.h_poca_xy.SetXTitle("x [cm]")
153 self.h_poca_xy.SetYTitle("y [cm]")
154
155 self.h_poca_z = TH1F("poca_z", "POCA distribution in z", 100, -2.0, 2.0)
156 self.h_poca_z.SetXTitle("z [cm]")
157
158 self.h_hitsCDC = TH1F("cdc_hits", "CDC hits", 100, 0.0, 100.0)
159 self.h_hitsCDC.SetXTitle("number of CDC hits")
160
161 self.h_Ecms = TH1F("Ecms", "c.m.s. energy of muon", 300, 4.5, 6.0)
162 self.h_Ecms.SetXTitle("E_{cm} [GeV/c]")
163
164
165 self.tree = ROOT.TTree('tree', 'time resolution')
166
167 self.data = TreeStruct()
168 self.data.itrk = 0
169
170 for key in TreeStruct.__dict__.keys():
171 if '__' not in key:
172 formstring = '/F'
173 if isinstance(self.data.__getattribute__(key), int):
174 formstring = '/I'
175 if key in int_arrays:
176 formstring = '[nfot]/I'
177 elif key in float_arrays:
178 formstring = '[nfot]/F'
179 self.tree.Branch(key, addressof(self.data, key), key + formstring)
180
Utility for the track selection - used in various calibration modules.
Definition: TrackSelector.h:27

◆ pdfHistogram()

def pdfHistogram (   self,
  pdf 
)
 make histogram of PDF peak positions for the first 20 tracks 

Definition at line 187 of file cdst_timeResoNtuple.py.

187 def pdfHistogram(self, pdf):
188 ''' make histogram of PDF peak positions for the first 20 tracks '''
189
190 self.nhisto += 1
191 if self.nhisto > 20:
192 return
193 h = TH2F('pdf' + str(self.nhisto), 'muon PDF, itrk = ' + str(self.data.itrk),
194 512, 0.0, 512.0, 1000, 0.0, 75.0)
195 h.SetXTitle('pixelID - 1')
196 h.SetYTitle('peak positions [ns]')
197 for x in range(512):
198 peaks = pdf[x]
199 for peak in peaks:
200 h.Fill(x, peak.mean)
201 h.Write()
202

◆ sortPeaks()

def sortPeaks (   self,
  unsortedPeaks 
)
 sort PDF peaks according to their positions 

Definition at line 181 of file cdst_timeResoNtuple.py.

181 def sortPeaks(self, unsortedPeaks):
182 ''' sort PDF peaks according to their positions '''
183
184 py_list = [x for x in unsortedPeaks]
185 return sorted(py_list, key=lambda x: (x.mean))
186

◆ terminate()

def terminate (   self)
 terminate: close root file 

Definition at line 351 of file cdst_timeResoNtuple.py.

351 def terminate(self):
352 ''' terminate: close root file '''
353
354 self.file.cd()
355 self.file.Write()
356 self.file.Close()
357
358
359# Create path

Member Data Documentation

◆ bunchOffset

bunchOffset

histogram of bunch offset

Definition at line 124 of file cdst_timeResoNtuple.py.

◆ data

data

data structure

Definition at line 167 of file cdst_timeResoNtuple.py.

◆ file

file

file object

Definition at line 121 of file cdst_timeResoNtuple.py.

◆ h_charge

h_charge

histogram of charge

Definition at line 147 of file cdst_timeResoNtuple.py.

◆ h_cth

h_cth

histogram of cos(theta)

Definition at line 135 of file cdst_timeResoNtuple.py.

◆ h_cth_vs_p

h_cth_vs_p

histogram cos(theta) vs, momentum

Definition at line 127 of file cdst_timeResoNtuple.py.

◆ h_Ecms

h_Ecms

histogram of Ecms

Definition at line 161 of file cdst_timeResoNtuple.py.

◆ h_hitsCDC

h_hitsCDC

histogram of number of hits in CDC

Definition at line 158 of file cdst_timeResoNtuple.py.

◆ h_momentum

h_momentum

histogram of momentum

Definition at line 132 of file cdst_timeResoNtuple.py.

◆ h_phi

h_phi

histogram of local phi

Definition at line 138 of file cdst_timeResoNtuple.py.

◆ h_poca_xy

h_poca_xy

histogram of POCA x-y

Definition at line 150 of file cdst_timeResoNtuple.py.

◆ h_poca_z

h_poca_z

histogram of POCA z

Definition at line 155 of file cdst_timeResoNtuple.py.

◆ h_x

h_x

histogram of local x

Definition at line 144 of file cdst_timeResoNtuple.py.

◆ h_z

h_z

histogram of local z

Definition at line 141 of file cdst_timeResoNtuple.py.

◆ nhisto

int nhisto = 0
static

histogram counter

Definition at line 103 of file cdst_timeResoNtuple.py.

◆ tree

tree

tree object

Definition at line 165 of file cdst_timeResoNtuple.py.

◆ trk_selector

trk_selector

track selector - selection of muons from ee->mumu events

Definition at line 115 of file cdst_timeResoNtuple.py.


The documentation for this class was generated from the following file: