Belle II Software  release-08-01-10
cdst_timeResoNtuple.py
1 #!/usr/bin/env python3
2 # -*- coding: utf-8 -*-
3 
4 
11 
12 # ---------------------------------------------------------------------------------------
13 # Makes ntuple w/ variable size arrays to study single photon time resolution with di-muon data
14 # by comparing photon times w.r.t leading PDF peak position.
15 # Provides also full information of a PDF peak assocoated with a detected photon.
16 # The following draw command must show a peak that represents an overall time resolution:
17 # tree->Draw("time-t0>>h(100,-1,2)", "t0<10")
18 #
19 # Usage: basf2 cdst_timeResoNtuple.py -i <cdst_file-dimuon_skim.root>
20 # ---------------------------------------------------------------------------------------
21 
22 import basf2 as b2
23 from ROOT import Belle2
24 from ROOT import TH1F, TH2F
25 from ROOT import addressof, gROOT
26 import math
27 import ROOT
28 
29 # global tags
30 # ******************************************************************************************************************
31 # note: The patching global tags and their order are bucket number and basf2 version dependent.
32 # Given below is what is needed for cdst files of bucket 16 calibration and January-2023 development version.
33 # ******************************************************************************************************************
34 b2.conditions.override_globaltags()
35 b2.conditions.append_globaltag('patch_main_release-07_noTOP')
36 b2.conditions.append_globaltag('data_reprocessing_proc13') # experiments 7 - 18
37 # b2.conditions.append_globaltag('data_reprocessing_prompt') # experiments 20 - 26
38 b2.conditions.append_globaltag('online')
39 
40 
41 gROOT.ProcessLine('struct TreeStruct {\
42  int run; /* run number */ \
43  int itrk; /* track counter */ \
44  float offset; /* current bunch offset */ \
45  int usedTrk; /* number of tracks used for bunch reconstruction */ \
46  int slot; /* slot ID */ \
47  float p; /* extHit momentum */ \
48  float cth; /* extHit local cos(theta) (global must be pretty much the same) */ \
49  float phi; /* extHit local phi */ \
50  float z; /* extHit local z */ \
51  float x; /* extHit local x */ \
52  float tof; /* extHit time-of-flight */ \
53  int charge; /* particle charge */ \
54  float poca_x; /* POCA x */ \
55  float poca_y; /* POCA y */ \
56  float poca_z; /* POCA z */ \
57  int hitsCDC; /* number of CDC hits */ \
58  float Ecms; /* CMS energy of muon */ \
59  float rec_t0; /* t0 of this track determined with bunch finder */ \
60  int valid_t0; /* reconstruction status of t0 */ \
61  int nfot; /* number of photons */ \
62  int channel[1000]; /* channel */ \
63  int pixel[1000] ; /* pixel ID */ \
64  float time[1000]; /* photon time */ \
65  float timeErr[1000]; /* photon time uncertainty */ \
66  int pulseHeight[1000]; /* photon pulse height */ \
67  float pulseWidth[1000]; /* photon pulse width */ \
68  int sample[1000]; /* sample number modulo 256 */ \
69  int status[1000]; /* calibration status bits */ \
70  float t0[1000]; /* leading PDF peak: position */ \
71  float wid0[1000]; /* leading PDF peak: width w/o TTS */ \
72  float t1[1000]; /* next to leading PDF peak: position */ \
73  float t_pdf[1000]; /* associated pdf peak time */ \
74  float wid[1000]; /* associated pdf peak width */ \
75  int type[1000]; /* 0 bkg, 1 direct, 2 reflected */ \
76  int nx[1000]; /* total number of reflections in x */ \
77  int ny[1000]; /* total number of reflections in y */ \
78  int nxm[1000]; /* number of reflections in x before mirror */ \
79  int nym[1000]; /* number of reflections in y before mirror */ \
80  int nxe[1000]; /* number of reflections in x in prism */ \
81  int nye[1000]; /* number of reflections in y in prism */ \
82  int nys[1000]; /* number of reflections on slanted surface of prism */ \
83  float fic[1000]; /* Cherenkov azimuthal angle in degrees */ \
84  float xd[1000]; /* unfolded x coordinate of a pixel */ \
85  float yd[1000]; /* unfolded y coordinate of a pixel */ \
86  float xm[1000]; /* unfolded x coordinate of a reconstructed point on mirror */ \
87  float kx[1000]; /* reconstructed photon direction in x at emission */ \
88  float ky[1000]; /* reconstructed photon direction in y at emission */ \
89  float alpha[1000]; /* impact angle on photo-cathode [degrees] */ \
90 };')
91 
92 int_arrays = ['channel', 'pixel', 'pulseHeight', 'sample', 'status', 'type', 'nx', 'ny',
93  'nxm', 'nym', 'nxe', 'nye', 'nys']
94 float_arrays = ['time', 'timeErr', 'pulseWidth', 't0', 'wid0', 't1', 't_pdf', 'wid', 'alpha', 'fic',
95  'xd', 'yd', 'xm', 'kx', 'ky']
96 
97 from ROOT import TreeStruct # noqa
98 
99 
100 class Ntuple(b2.Module):
101  ''' Makes a flat ntuple '''
102 
103 
104  nhisto = 0
105 
106  def initialize(self):
107  ''' initialize: open root file, construct ntuple '''
108 
109  evtMetaData = Belle2.PyStoreObj('EventMetaData')
110  expNo = evtMetaData.obj().getExperiment()
111  runNo = evtMetaData.obj().getRun()
112  exp_run = '-e' + '{:0=4d}'.format(expNo) + '-r' + '{:0=5d}'.format(runNo)
113  outName = 'timeResoNtuple' + exp_run + '.root'
114 
115 
116  self.trk_selectortrk_selector = Belle2.TOP.TrackSelector("dimuon")
117  self.trk_selectortrk_selector.setDeltaEcms(0.1)
118  self.trk_selectortrk_selector.setCutOnPOCA(2.0, 4.0)
119  self.trk_selectortrk_selector.setCutOnLocalZ(-130.0, 130.0)
120 
121 
122  self.filefile = ROOT.TFile(outName, 'recreate')
123 
124 
125  self.bunchOffsetbunchOffset = TH1F("bunchOffset", "bunch offset", 100, -1.0, 1.0)
126  self.bunchOffsetbunchOffset.SetXTitle("bunch offset [ns]")
127 
128  self.h_cth_vs_ph_cth_vs_p = TH2F("cth_vs_p", "local cos #theta vs. p", 100, 0.0, 10.0,
129  100, -1.0, 1.0)
130  self.h_cth_vs_ph_cth_vs_p.SetXTitle("p [GeV/c]")
131  self.h_cth_vs_ph_cth_vs_p.SetYTitle("cos #theta")
132 
133  self.h_momentumh_momentum = TH1F("momentum", "momentum", 100, 0.0, 10.0)
134  self.h_momentumh_momentum.SetXTitle("p [GeV/c]")
135 
136  self.h_cthh_cth = TH1F("cos_theta", "local cos #theta", 100, -1.0, 1.0)
137  self.h_cthh_cth.SetXTitle("cos #theta")
138 
139  self.h_phih_phi = TH1F("local_phi", "local phi", 100, -math.pi, math.pi)
140  self.h_phih_phi.SetXTitle("local #phi")
141 
142  self.h_zh_z = TH1F("local_z", "local z", 100, -140.0, 140.0)
143  self.h_zh_z.SetXTitle("local z [cm]")
144 
145  self.h_xh_x = TH1F("local_x", "local x", 100, -23.0, 23.0)
146  self.h_xh_x.SetXTitle("local x [cm]")
147 
148  self.h_chargeh_charge = TH1F("charge", "charge", 3, -1.5, 1.5)
149  self.h_chargeh_charge.SetXTitle("charge")
150 
151  self.h_poca_xyh_poca_xy = TH2F("poca_xy", "POCA distribution in x-y", 100, -1.0, 1.0,
152  100, -1.0, 1.0)
153  self.h_poca_xyh_poca_xy.SetXTitle("x [cm]")
154  self.h_poca_xyh_poca_xy.SetYTitle("y [cm]")
155 
156  self.h_poca_zh_poca_z = TH1F("poca_z", "POCA distribution in z", 100, -2.0, 2.0)
157  self.h_poca_zh_poca_z.SetXTitle("z [cm]")
158 
159  self.h_hitsCDCh_hitsCDC = TH1F("cdc_hits", "CDC hits", 100, 0.0, 100.0)
160  self.h_hitsCDCh_hitsCDC.SetXTitle("number of CDC hits")
161 
162  self.h_Ecmsh_Ecms = TH1F("Ecms", "c.m.s. energy of muon", 300, 4.5, 6.0)
163  self.h_Ecmsh_Ecms.SetXTitle("E_{cm} [GeV/c]")
164 
165 
166  self.treetree = ROOT.TTree('tree', 'time resolution')
167 
168  self.datadata = TreeStruct()
169  self.datadata.itrk = 0
170 
171  for key in TreeStruct.__dict__.keys():
172  if '__' not in key:
173  formstring = '/F'
174  if isinstance(self.datadata.__getattribute__(key), int):
175  formstring = '/I'
176  if key in int_arrays:
177  formstring = '[nfot]/I'
178  elif key in float_arrays:
179  formstring = '[nfot]/F'
180  self.treetree.Branch(key, addressof(self.datadata, key), key + formstring)
181 
182  def sortPeaks(self, unsortedPeaks):
183  ''' sort PDF peaks according to their positions '''
184 
185  py_list = [x for x in unsortedPeaks]
186  return sorted(py_list, key=lambda x: (x.mean))
187 
188  def pdfHistogram(self, pdf):
189  ''' make histogram of PDF peak positions for the first 20 tracks '''
190 
191  self.nhistonhisto += 1
192  if self.nhistonhisto > 20:
193  return
194  h = TH2F('pdf' + str(self.nhistonhisto), 'muon PDF, itrk = ' + str(self.datadata.itrk),
195  512, 0.0, 512.0, 1000, 0.0, 75.0)
196  h.SetXTitle('pixelID - 1')
197  h.SetYTitle('peak positions [ns]')
198  for x in range(512):
199  peaks = pdf[x]
200  for peak in peaks:
201  h.Fill(x, peak.mean)
202  h.Write()
203 
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.bunchOffsetbunchOffset.Fill(recBunch.getCurrentOffset())
214 
215  evtMetaData = Belle2.PyStoreObj('EventMetaData')
216  self.datadata.run = evtMetaData.getRun()
217  self.datadata.offset = recBunch.getCurrentOffset()
218  self.datadata.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_selectortrk_selector.isSelected(trk):
225  continue
226  pdfs = track.getRelated('TOPPDFCollections')
227  if not pdfs:
228  continue
229  self.datadata.slot = pdfs.getModuleID()
230  momentum = pdfs.getAssociatedLocalMomentum()
231  position = pdfs.getAssociatedLocalHit()
232  self.datadata.p = momentum.R()
233  self.datadata.cth = math.cos(momentum.Theta())
234  self.datadata.phi = momentum.Phi()
235  self.datadata.z = position.Z()
236  self.datadata.x = position.X()
237  self.datadata.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.datadata.charge = tfit.getChargeSign()
244  pocaPosition = tfit.getPosition()
245  self.datadata.poca_x = pocaPosition.X()
246  self.datadata.poca_y = pocaPosition.Y()
247  self.datadata.poca_z = pocaPosition.Z()
248  self.datadata.hitsCDC = tfit.getHitPatternCDC().getNHits()
249  self.datadata.Ecms = self.trk_selectortrk_selector.getCMSEnergy()
250  try:
251  extHit = trk.getExtHit()
252  timeZero = extHit.getRelated('TOPTimeZeros')
253  self.datadata.rec_t0 = timeZero.getTime()
254  self.datadata.valid_t0 = timeZero.isValid()
255  except BaseException:
256  self.datadata.rec_t0 = 0
257  self.datadata.valid_t0 = 0
258 
259  self.h_cth_vs_ph_cth_vs_p.Fill(self.datadata.p, self.datadata.cth)
260  self.h_momentumh_momentum.Fill(self.datadata.p)
261  self.h_cthh_cth.Fill(self.datadata.cth)
262  self.h_phih_phi.Fill(self.datadata.phi)
263  self.h_zh_z.Fill(self.datadata.z)
264  self.h_xh_x.Fill(self.datadata.x)
265  self.h_chargeh_charge.Fill(self.datadata.charge)
266  self.h_poca_xyh_poca_xy.Fill(self.datadata.poca_x, self.datadata.poca_y)
267  self.h_poca_zh_poca_z.Fill(self.datadata.poca_z)
268  self.h_hitsCDCh_hitsCDC.Fill(self.datadata.hitsCDC)
269  self.h_Ecmsh_Ecms.Fill(self.datadata.Ecms)
270  try:
271  pdf = pdfs.getHypothesisPDF(13)
272  except BaseException:
273  b2.B2ERROR("No PDF available for PDG = 13")
274  continue
275  self.datadata.itrk += 1
276  self.pdfHistogrampdfHistogram(pdf)
277  emi = trk.getEmissionPoint()
278  x0 = emi.position.X()
279  z0 = emi.position.Y()
280  self.datadata.nfot = 0
281  for digit in Belle2.PyStoreArray('TOPDigits'):
282  if digit.getModuleID() == self.datadata.slot and digit.getHitQuality() == 1:
283  k = self.datadata.nfot
284  if k >= 1000:
285  continue
286  self.datadata.nfot += 1
287  self.datadata.channel[k] = digit.getChannel()
288  self.datadata.pixel[k] = digit.getPixelID()
289  self.datadata.time[k] = digit.getTime()
290  self.datadata.timeErr[k] = digit.getTimeError()
291  self.datadata.pulseHeight[k] = digit.getPulseHeight()
292  self.datadata.pulseWidth[k] = digit.getPulseWidth()
293  self.datadata.sample[k] = digit.getModulo256Sample()
294  self.datadata.status[k] = digit.getStatus()
295  peaks = pdf[digit.getPixelID() - 1]
296  if peaks.empty():
297  self.datadata.t0[k] = 0
298  self.datadata.wid0[k] = 0
299  self.datadata.t1[k] = 0
300  else:
301  sorted_peaks = self.sortPeakssortPeaks(peaks)
302  self.datadata.t0[k] = sorted_peaks[0].mean
303  self.datadata.wid0[k] = sorted_peaks[0].width
304  self.datadata.t1[k] = self.datadata.t0[k] + 100
305  if peaks.size() > 1:
306  self.datadata.t1[k] = sorted_peaks[1].mean
307  self.datadata.t_pdf[k] = 0
308  self.datadata.wid[k] = 0
309  self.datadata.fic[k] = 0
310  self.datadata.type[k] = 0
311  self.datadata.nx[k] = 0
312  self.datadata.ny[k] = 0
313  self.datadata.nxm[k] = 0
314  self.datadata.nym[k] = 0
315  self.datadata.nxe[k] = 0
316  self.datadata.nye[k] = 0
317  self.datadata.nys[k] = 0
318  self.datadata.xd[k] = 0
319  self.datadata.yd[k] = 0
320  self.datadata.xm[k] = 0
321  self.datadata.kx[k] = 0
322  self.datadata.ky[k] = 0
323  self.datadata.alpha[k] = 0
324  assocPDF = digit.getRelated('TOPAssociatedPDFs')
325  if assocPDF:
326  pik = assocPDF.getSinglePeak()
327  if pik:
328  self.datadata.t_pdf[k] = pik.position
329  self.datadata.wid[k] = pik.width
330  self.datadata.fic[k] = math.degrees(pik.fic)
331  self.datadata.type[k] = pik.type
332  self.datadata.nx[k] = pik.nx
333  self.datadata.ny[k] = pik.ny
334  self.datadata.nxm[k] = pik.nxm
335  self.datadata.nym[k] = pik.nym
336  self.datadata.nxe[k] = pik.nxe
337  self.datadata.nye[k] = pik.nye
338  if pik.kyd > 0:
339  self.datadata.nys[k] = int(pik.nye / 2)
340  else:
341  self.datadata.nys[k] = int((pik.nye + 1) / 2)
342  self.datadata.xd[k] = pik.xd
343  self.datadata.yd[k] = pik.yd
344  if pik.kze > 0:
345  self.datadata.xm[k] = x0 + pik.kxe / pik.kze * (130.0 - z0)
346  self.datadata.kx[k] = pik.kxe
347  self.datadata.ky[k] = pik.kye
348  self.datadata.alpha[k] = math.degrees(math.acos(abs(pik.kzd)))
349 
350  self.treetree.Fill()
351 
352  def terminate(self):
353  ''' terminate: close root file '''
354 
355  self.filefile.cd()
356  self.filefile.Write()
357  self.filefile.Close()
358 
359 
360 # Create path
361 main = b2.create_path()
362 
363 # Input: cdst file(s), use -i option
364 main.add_module('RootInput')
365 
366 # Reconstruction: TOP only
367 main.add_module('Gearbox')
368 main.add_module('Geometry')
369 main.add_module('Ext')
370 main.add_module('TOPUnpacker')
371 main.add_module('TOPRawDigitConverter')
372 main.add_module('TOPChannelMasker')
373 main.add_module('TOPBunchFinder')
374 
375 # Make a muon PDF available at datastore
376 main.add_module('TOPPDFDebugger', pdgCodes=[13]) # default
377 
378 # Write ntuple
379 main.add_module(Ntuple())
380 
381 # Print progress
382 progress = b2.register_module('Progress')
383 main.add_module(progress)
384 
385 # Process events
386 b2.process(main)
387 
388 # Print call statistics
389 print(b2.statistics)
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
Utility for the track selection - used in various calibration modules.
Definition: TrackSelector.h:27
trk_selector
track selector - selection of muons from ee->mumu events
def sortPeaks(self, unsortedPeaks)
h_cth
histogram of cos(theta)
h_cth_vs_p
histogram cos(theta) vs, momentum
h_phi
histogram of local phi
h_hitsCDC
histogram of number of hits in CDC
h_poca_xy
histogram of POCA x-y
h_momentum
histogram of momentum
bunchOffset
histogram of bunch offset
static ExpRun getRun(std::map< ExpRun, std::pair< double, double >> runs, double t)
Get exp number + run number from time.
Definition: Splitter.cc:262