Belle II Software  release-08-01-10
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 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 100 of file cdst_timeResoNtuple.py.

Member Function Documentation

◆ event()

def event (   self)
 event processing 

Definition at line 204 of file cdst_timeResoNtuple.py.

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

◆ pdfHistogram()

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

Definition at line 188 of file cdst_timeResoNtuple.py.

◆ sortPeaks()

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

Definition at line 182 of file cdst_timeResoNtuple.py.

◆ terminate()

def terminate (   self)
 terminate: close root file 

Definition at line 352 of file cdst_timeResoNtuple.py.


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