Belle II Software development
plot_LL_diff.py
1#!/usr/bin/env python3
2
3
10
11
20
21import basf2 as b2
22
23from ROOT import Belle2
24from ROOT import TCanvas, TH2F
25from simulation import add_simulation
26from reconstruction import add_reconstruction
27
28# the histograms to fill
29hist = [TH2F(
30 'hist' + str(i),
31 'hist' + str(i),
32 50,
33 0.0,
34 3.5,
35 50,
36 -100,
37 100,
38) for i in range(4)]
39
40hist[0].SetTitle('True pions, dE/dx')
41hist[1].SetTitle('True kaons, dE/dx')
42hist[2].SetTitle('True pions, TOP')
43hist[3].SetTitle('True kaons, TOP')
44
45# define the python module to extract and plot PID information
46
47
48class MinModule(b2.Module):
49
50 """
51 Get LL differences from PIDLikelihood and fill them into histograms
52 """
53
54 def __init__(self):
55 """
56 call constructor of base class, required.
57 """
58
59 super().__init__()
60
61 def event(self):
62 """
63 reimplement Module::event()
64 """
65
66 pids = Belle2.PyStoreArray('PIDLikelihoods')
67 for pid in pids:
68 track = pid.getRelatedFrom('Tracks')
69 if track:
70 mcpart = track.getRelatedTo('MCParticles')
71 if not track or not mcpart:
72 # some tracks don't have an mcparticle
73 b2.B2WARNING('problems with track <-> mcparticle relations')
74 event = Belle2.PyStoreObj('EventMetaData').obj().getEvent()
75 print(f'event: {int(event)}, track: {int(track.getArrayIndex())}')
76 else:
77 pdg = abs(mcpart.getPDG())
78 momentumVec = mcpart.getMomentum()
79 momentum = momentumVec.Mag()
80 if momentum > 3.5: # cut off
81 continue
82 # theta = momentumVec.CosTheta()
83
84 # particle to compare with pions
85 selectedpart = Belle2.Const.kaon
86 pid_dedx = Belle2.Const.PIDDetectorSet(Belle2.Const.CDC)
87 pid_top = Belle2.Const.PIDDetectorSet(Belle2.Const.TOP)
88 logl_sel = pid.getLogL(selectedpart, pid_dedx)
89 logl_pi = pid.getLogL(Belle2.Const.pion, pid_dedx)
90 dedx_DLL = logl_pi - logl_sel
91
92 logl_sel = pid.getLogL(selectedpart, pid_top)
93 logl_pi = pid.getLogL(Belle2.Const.pion, pid_top)
94 top_DLL = logl_pi - logl_sel
95
96 if pdg == selectedpart.getPDGCode():
97 hist[0].Fill(momentum, dedx_DLL)
98 hist[2].Fill(momentum, top_DLL)
99 elif pdg == 211:
100 hist[1].Fill(momentum, dedx_DLL)
101 hist[3].Fill(momentum, top_DLL)
102
103 def terminate(self):
104 """
105 Draw histograms on canvas and save image.
106 """
107
108 c1 = TCanvas('c', 'c', 800, 600)
109 c1.Divide(2, 2)
110 c1.Show()
111 for i in range(4):
112 c1.cd(i + 1)
113 hist[i].Draw('colz')
114 # c1.Draw()
115
116 c1.SaveAs('ll_diff.png')
117 print('Output plots saved in ll_diff.png')
118
119
120# create path
121main = b2.create_path()
122
123# define number of events
124eventinfosetter = b2.register_module('EventInfoSetter')
125eventinfosetter.param('evtNumList', [20])
126main.add_module(eventinfosetter)
127main.add_module('EvtGenInput')
128
129# do full simulation and reconstruction
130add_simulation(main)
131add_reconstruction(main)
132
133# draw histograms of log likelihood differences
134main.add_module(MinModule())
135
136# process events and print call statistics
137b2.process(main)
138print(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