Belle II Software development
PhokharaEvtgenAnalyze.py
1#!/usr/bin/env python3
2
3
10
11"""
12<header>
13 <input>PhokharaEvtgenData.root</input>
14 <output>PhokharaEvtgenAnalysis.root</output>
15 <contact>Kirill Chilikin (K.A.Chilikin@inp.nsk.su)</contact>
16 <description>Analysis of e+ e- -> J/psi eta_c events.</description>
17</header>
18"""
19
20import basf2 as b2
21import ROOT
22import numpy
23from ROOT import Belle2
24
25
27 """ Analysis module for PhokharaEvtgen. """
28
29 def __init__(self):
30 """Initialization."""
31 super().__init__()
32
33 self.output_file = ROOT.TFile('PhokharaEvtgenAnalysis.root', 'recreate')
34
35 self.tree = ROOT.TTree('tree', '')
36
37 self.ecms = numpy.zeros(1, dtype=numpy.float32)
38
39 self.gamma_e = numpy.zeros(1, dtype=numpy.float32)
40
41 self.gamma_px = numpy.zeros(1, dtype=numpy.float32)
42
43 self.gamma_py = numpy.zeros(1, dtype=numpy.float32)
44
45 self.gamma_pz = numpy.zeros(1, dtype=numpy.float32)
46
47 self.jpsi_e = numpy.zeros(1, dtype=numpy.float32)
48
49 self.jpsi_px = numpy.zeros(1, dtype=numpy.float32)
50
51 self.jpsi_py = numpy.zeros(1, dtype=numpy.float32)
52
53 self.jpsi_pz = numpy.zeros(1, dtype=numpy.float32)
54
55 self.lepton_e = numpy.zeros(1, dtype=numpy.float32)
56
57 self.lepton_px = numpy.zeros(1, dtype=numpy.float32)
58
59 self.lepton_py = numpy.zeros(1, dtype=numpy.float32)
60
61 self.lepton_pz = numpy.zeros(1, dtype=numpy.float32)
62 self.tree.Branch('ecms', self.ecms, 'ecms/F')
63 self.tree.Branch('gamma_e', self.gamma_e, 'gamma_e/F')
64 self.tree.Branch('gamma_px', self.gamma_px, 'gamma_px/F')
65 self.tree.Branch('gamma_py', self.gamma_py, 'gamma_py/F')
66 self.tree.Branch('gamma_pz', self.gamma_pz, 'gamma_pz/F')
67 self.tree.Branch('jpsi_e', self.jpsi_e, 'jpsi_e/F')
68 self.tree.Branch('jpsi_px', self.jpsi_px, 'jpsi_px/F')
69 self.tree.Branch('jpsi_py', self.jpsi_py, 'jpsi_py/F')
70 self.tree.Branch('jpsi_pz', self.jpsi_pz, 'jpsi_pz/F')
71 self.tree.Branch('lepton_e', self.lepton_e, 'lepton_e/F')
72 self.tree.Branch('lepton_px', self.lepton_px, 'lepton_px/F')
73 self.tree.Branch('lepton_py', self.lepton_py, 'lepton_py/F')
74 self.tree.Branch('lepton_pz', self.lepton_pz, 'lepton_pz/F')
75
76 def event(self):
77 """ Event function. """
78 mc_initial_particles = Belle2.PyStoreObj('MCInitialParticles')
79 self.ecms[0] = mc_initial_particles.getMass()
80 mc_particles = Belle2.PyStoreArray('MCParticles')
81 for mc_particle in mc_particles:
82 # Select virtual photons.
83 if (mc_particle.getPDG() != 10022):
84 continue
85 p = mc_particle.getMomentum()
86 self.gamma_e[0] = mc_particle.getEnergy()
87 self.gamma_px[0] = p.X()
88 self.gamma_py[0] = p.Y()
89 self.gamma_pz[0] = p.Z()
90 jpsi = mc_particle.getDaughters()[0]
91 p = jpsi.getMomentum()
92 self.jpsi_e[0] = jpsi.getEnergy()
93 self.jpsi_px[0] = p.X()
94 self.jpsi_py[0] = p.Y()
95 self.jpsi_pz[0] = p.Z()
96 lepton = jpsi.getDaughters()[0]
97 p = lepton.getMomentum()
98 self.lepton_e[0] = lepton.getEnergy()
99 self.lepton_px[0] = p.X()
100 self.lepton_py[0] = p.Y()
101 self.lepton_pz[0] = p.Z()
102 self.tree.Fill()
103
104 def terminate(self):
105 """ Termination function. """
106 self.output_file.cd()
107 self.tree.Write()
108 self.output_file.Close()
109
110
111# Input.
112root_input = b2.register_module('RootInput')
113root_input.param('inputFileName', 'PhokharaEvtgenData.root')
114
115# Analysis.
116phokhara_evtgen = PhokharaEvtgenAnalysisModule()
117
118# Create main path.
119main = b2.create_path()
120
121# Add modules to main path
122main.add_module(root_input)
123main.add_module(phokhara_evtgen)
124
125main.add_module('Progress')
126# Run.
127b2.process(main)
A (simplified) python wrapper for StoreArray.
a (simplified) python wrapper for StoreObjPtr.
Definition PyStoreObj.h:67