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