Belle II Software  release-05-02-19
Ntuple Class Reference
Inheritance diagram for Ntuple:
Collaboration diagram for Ntuple:

Public Member Functions

def initialize (self)
 
def sortPeaks (self, unsortedPeaks)
 
def pdfHistogram (self, pdf)
 
def getEcms (self, tfit, chargedStable)
 
def trackSelected (self)
 
def getTOF (self, track, pdg, slot)
 
def event (self)
 
def terminate (self)
 

Public Attributes

 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 87 of file cdst_timeResoNtuple.py.

Member Function Documentation

◆ event()

def event (   self)
event processing 

Definition at line 232 of file cdst_timeResoNtuple.py.

232  def event(self):
233  ''' event processing '''
234 
235  recBunch = Belle2.PyStoreObj('TOPRecBunch')
236  if not recBunch:
237  B2ERROR('no TOPRecBunch')
238  return
239  if not recBunch.isReconstructed():
240  return
241  self.bunchOffset.Fill(recBunch.getCurrentOffset())
242 
243  evtMetaData = Belle2.PyStoreObj('EventMetaData')
244  self.data.run = evtMetaData.getRun()
245  self.data.offset = recBunch.getCurrentOffset()
246  self.data.usedTrk = recBunch.getUsedTracks()
247 
248  for track in Belle2.PyStoreArray('Tracks'):
249  pdfs = track.getRelated('TOPPDFCollections')
250  if not pdfs:
251  continue
252  self.data.slot = pdfs.getModuleID()
253  momentum = pdfs.getAssociatedLocalMomentum()
254  position = pdfs.getAssociatedLocalHit()
255  self.data.p = momentum.Mag()
256  self.data.cth = momentum.CosTheta()
257  self.data.phi = momentum.Phi()
258  self.data.z = position.Z()
259  self.data.x = position.X()
260  self.data.tof = self.getTOF(track, 13, self.data.slot)
261  try:
262  tfit = track.getTrackFitResultWithClosestMass(Belle2.Const.muon)
263  except:
264  B2ERROR("No trackFitResult available")
265  continue
266  self.data.charge = tfit.getChargeSign()
267  pocaPosition = tfit.getPosition()
268  self.data.poca_x = pocaPosition.X()
269  self.data.poca_y = pocaPosition.Y()
270  self.data.poca_z = pocaPosition.Z()
271  self.data.hitsCDC = tfit.getHitPatternCDC().getNHits()
272  Ecms = self.getEcms(tfit, Belle2.Const.muon)
273  self.data.dEcms = Ecms[0] - Ecms[1]
274  if not self.trackSelected():
275  continue
276  try:
277  topll = track.getRelated('TOPLikelihoods')
278  extHit = topll.getRelated('ExtHits')
279  timeZero = extHit.getRelated('TOPTimeZeros')
280  self.data.rec_t0 = timeZero.getTime()
281  self.data.valid_t0 = timeZero.isValid()
282  except:
283  self.data.rec_t0 = 0
284  self.data.valid_t0 = 0
285 
286  self.h_cth_vs_p.Fill(self.data.p, self.data.cth)
287  self.h_momentum.Fill(self.data.p)
288  self.h_cth.Fill(self.data.cth)
289  self.h_phi.Fill(self.data.phi)
290  self.h_z.Fill(self.data.z)
291  self.h_x.Fill(self.data.x)
292  self.h_charge.Fill(self.data.charge)
293  self.h_poca_xy.Fill(self.data.poca_x, self.data.poca_y)
294  self.h_poca_z.Fill(self.data.poca_z)
295  self.h_hitsCDC.Fill(self.data.hitsCDC)
296  self.h_Ecms.Fill(Ecms[0])
297  try:
298  pdf = pdfs.getHypothesisPDF(13)
299  except:
300  B2ERROR("No PDF available for PDG = 13")
301  continue
302  self.data.itrk += 1
303  self.pdfHistogram(pdf)
304  x0 = position.X() - momentum.X() / momentum.Y() * position.Y() # emission
305  z0 = position.Z() - momentum.Z() / momentum.Y() * position.Y() # emission
306  self.data.nfot = 0
307  for digit in Belle2.PyStoreArray('TOPDigits'):
308  if digit.getModuleID() == self.data.slot and digit.getHitQuality() == 1:
309  k = self.data.nfot
310  if k >= 1000:
311  continue
312  self.data.nfot += 1
313  self.data.channel[k] = digit.getChannel()
314  self.data.pixel[k] = digit.getPixelID()
315  self.data.time[k] = digit.getTime()
316  self.data.timeErr[k] = digit.getTimeError()
317  self.data.pulseHeight[k] = digit.getPulseHeight()
318  self.data.pulseWidth[k] = digit.getPulseWidth()
319  self.data.sample[k] = digit.getModulo256Sample()
320  self.data.status[k] = digit.getStatus()
321  peaks = pdf[digit.getPixelID() - 1]
322  if peaks.empty():
323  self.data.t0[k] = 0
324  self.data.wid0[k] = 0
325  self.data.t1[k] = 0
326  else:
327  sorted_peaks = self.sortPeaks(peaks)
328  self.data.t0[k] = sorted_peaks[0].mean
329  self.data.wid0[k] = sorted_peaks[0].width
330  self.data.t1[k] = self.data.t0[k] + 100
331  if peaks.size() > 1:
332  self.data.t1[k] = sorted_peaks[1].mean
333  self.data.t_pdf[k] = 0
334  self.data.type[k] = 0
335  self.data.nx[k] = 0
336  self.data.ny[k] = 0
337  self.data.nxm[k] = 0
338  self.data.nym[k] = 0
339  self.data.nxe[k] = 0
340  self.data.nye[k] = 0
341  self.data.nys[k] = 0
342  self.data.xd[k] = 0
343  self.data.yd[k] = 0
344  self.data.xm[k] = 0
345  self.data.kx[k] = 0
346  self.data.ky[k] = 0
347  self.data.alpha[k] = 0
348  assocPDF = digit.getRelated('TOPAssociatedPDFs')
349  if assocPDF:
350  pik = assocPDF.getSinglePeak()
351  if pik:
352  self.data.t_pdf[k] = pik.position
353  self.data.type[k] = pik.type
354  self.data.nx[k] = pik.nx
355  self.data.ny[k] = pik.ny
356  self.data.nxm[k] = pik.nxm
357  self.data.nym[k] = pik.nym
358  self.data.nxe[k] = pik.nxe
359  self.data.nye[k] = pik.nye
360  if pik.kyd > 0:
361  self.data.nys[k] = int(pik.nye / 2)
362  else:
363  self.data.nys[k] = int((pik.nye + 1) / 2)
364  self.data.xd[k] = pik.xd
365  self.data.yd[k] = pik.yd
366  if pik.kze > 0:
367  self.data.xm[k] = x0 + pik.kxe / pik.kze * (130.0 - z0)
368  self.data.kx[k] = pik.kxe
369  self.data.ky[k] = pik.kye
370  self.data.alpha[k] = math.degrees(math.acos(abs(pik.kzd)))
371 
372  self.tree.Fill()
373 

◆ getEcms()

def getEcms (   self,
  tfit,
  chargedStable 
)
returns c.m.s. energy of a particle 

Definition at line 188 of file cdst_timeResoNtuple.py.

◆ getTOF()

def getTOF (   self,
  track,
  pdg,
  slot 
)
returns timo-of-flight of extrapolated track 

Definition at line 211 of file cdst_timeResoNtuple.py.

◆ initialize()

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

Definition at line 93 of file cdst_timeResoNtuple.py.

◆ pdfHistogram()

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

Definition at line 172 of file cdst_timeResoNtuple.py.

◆ sortPeaks()

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

Definition at line 166 of file cdst_timeResoNtuple.py.

◆ terminate()

def terminate (   self)
terminate: close root file 

Definition at line 374 of file cdst_timeResoNtuple.py.

◆ trackSelected()

def trackSelected (   self)
returns true if track fulfills selection criteria  

Definition at line 198 of file cdst_timeResoNtuple.py.


The documentation for this class was generated from the following file:
Belle2::PyStoreObj
a (simplified) python wrapper for StoreObjPtr.
Definition: PyStoreObj.h:69
ClusterEfficiency.ClusterEfficiency.event
def event(self)
Definition: ClusterEfficiency.py:146
Belle2::PyStoreArray
a (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:58