Belle II Software  release-08-01-10
pid_ttree_writer.py
1 #!/usr/bin/env python3
2 # -*- coding: utf-8 -*-
3 
4 
11 
12 
25 
26 import basf2 as b2
27 
28 # Some ROOT tools
29 import ROOT
30 from ROOT import Belle2
31 from ROOT import gROOT, AddressOf
32 
33 # Define a ROOT struct to hold output data in the TTree.
34 gROOT.ProcessLine('struct TreeStruct {\
35  float lld_dedx;\
36  float lld_top;\
37  float p;\
38  float phi;\
39  float costheta;\
40  float trackmomentum;\
41  int iskaon;\
42 };'
43  )
44 
45 
46 from ROOT import TreeStruct # noqa
47 # define the python module to save the PID information
48 
49 
50 class TreeWriterModule(b2.Module):
51 
52  """
53  This module writes its output to a ROOT tree.
54  Adapted from pxd/validation/PXDValidationTTreeSimHit.py
55  """
56 
57  def __init__(self):
58  """Initialize the module"""
59 
60  super(TreeWriterModule, self).__init__()
61 
62 
63  self.filefile = ROOT.TFile('PID_TTree.root', 'recreate')
64 
65  self.treetree = ROOT.TTree('tree', '')
66 
67  self.datadata = TreeStruct()
68  """ Declare tree branches """
69  for key in TreeStruct.__dict__.keys():
70  if '__' not in key:
71  formstring = '/F'
72  if isinstance(self.datadata.__getattribute__(key), int):
73  formstring = '/I'
74  self.treetree.Branch(key, AddressOf(self.datadata, key), key + formstring)
75 
76  def event(self):
77  """Store TOP and dE/dx info in tree"""
78 
79  pids = Belle2.PyStoreArray('PIDLikelihoods')
80  for pid in pids:
81  track = pid.getRelatedFrom('Tracks')
82  if track:
83  mcpart = track.getRelatedTo('MCParticles')
84  try:
85  pdg = abs(mcpart.getPDG())
86  if pdg != 211 and pdg != 321:
87  continue
88  momentumVec = mcpart.getMomentum()
89  momentum = momentumVec.Mag()
90  if momentum > 3.5: # cut off
91  continue
92  costheta = momentumVec.CosTheta()
93  phi = momentumVec.Phi()
94  fitresult = track.getTrackFitResultWithClosestMass(Belle2.Const.pion)
95  if fitresult:
96  trackmomentum = fitresult.getMomentum().Mag()
97  else:
98  trackmomentum = 0.0
99 
100  # particle to compare with pions
101  selectedpart = Belle2.Const.kaon
102  pid_dedx = Belle2.Const.PIDDetectorSet(Belle2.Const.CDC)
103  pid_top = Belle2.Const.PIDDetectorSet(Belle2.Const.TOP)
104  logl_sel = pid.getLogL(selectedpart, pid_dedx)
105  logl_pi = pid.getLogL(Belle2.Const.pion, pid_dedx)
106  dedx_DLL = logl_pi - logl_sel
107 
108  logl_sel = pid.getLogL(selectedpart, pid_top)
109  logl_pi = pid.getLogL(Belle2.Const.pion, pid_top)
110  top_DLL = logl_pi - logl_sel
111 
112  self.datadata.lld_dedx = dedx_DLL
113  self.datadata.lld_top = top_DLL
114  self.datadata.p = momentum
115  self.datadata.phi = phi
116  self.datadata.costheta = costheta
117  self.datadata.trackmomentum = trackmomentum
118  self.datadata.iskaon = pdg == 321
119 
120  # Fill tree
121  self.filefile.cd()
122  self.treetree.Fill()
123  except BaseException:
124 
125  # some tracks don't have an mcparticle (fixed now)
126  b2.B2WARNING('problems with track <-> mcparticle relations')
127  event = Belle2.PyStoreObj('EventMetaData').obj().getEvent()
128  print('event: %d, track: %d' % (event, track.getArrayIndex()))
129 
130  def terminate(self):
131  """ Close the output file."""
132 
133  self.filefile.cd()
134  self.filefile.Write()
135  self.filefile.Close()
136 
137 
138 # create path
139 main = b2.create_path()
140 
141 # use the input file defined via command line
142 main.add_module(b2.register_module('RootInput'))
143 
144 # add the python module defined above
145 main.add_module(TreeWriterModule())
146 
147 # process events and print call statistics
148 b2.process(main)
149 print(b2.statistics)
A class for sets of detector IDs whose content is limited to restricted set of valid detector IDs.
Definition: Const.h:287
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72
a (simplified) python wrapper for StoreObjPtr.
Definition: PyStoreObj.h:67
data
Instance of EventData class.