Belle II Software development
ParticleGunPlots.py
1#!/usr/bin/env python3
2
3
10
11import sys
12import math
13from basf2 import create_path, register_module, process, logging, \
14 LogLevel, Module, statistics
15
16# Load the required libraries
17from ROOT import TH1D, TH2D, TCanvas, Belle2, PyConfig
18
19logging.log_level = LogLevel.WARNING
20# reenable GUI thread for our canvas
21PyConfig.StartGuiThread = True
22
23# Create some histograms to be filled
24h_nTracks = TH1D('nTracks', 'Number of Tracks per Event', 100, 0, 50)
25h_pdg = TH1D('pid', 'PDG codes', 500, -250, 250)
26h_momentum = TH1D('momentum', 'Momentum p; p/GeV', 200, 0, 10)
27h_pt = TH1D('pt', 'Transverse Momentum p_{t};p_{t}/GeV', 200, 0, 10)
28h_phi = TH1D('phi', 'azimuthal angle #phi; #phi in degree', 200, -180, 180)
29h_theta = TH1D('theta', 'polar angle #theta; #theta in degree', 200, 0, 180)
30h_costheta = TH1D('costheta', 'cos(#theta); cos(#theta)', 200, -1, 1)
31h_xyvertex = TH2D('xyvertex', 'vertex in xy;x/#mum;y/#mum',
32 500, -500, 500, 500, -500, 500)
33h_zxvertex = TH2D('zxvertex', 'vertex in zx;z/#mum;x/#mum',
34 500, -500, 500, 500, -500, 500)
35h_zyvertex = TH2D('zyvertex', 'vertex in zy;z/#mum;y/#mum',
36 500, -500, 500, 500, -500, 500)
37
38
39class ShowMCParticles(Module):
40 """Simple module to collect some information about MCParticles"""
41
42 def event(self):
43 """Fill the histograms with the values of the MCParticle collection"""
44 mcParticles = Belle2.PyStoreArray('MCParticles')
45 h_nTracks.Fill(mcParticles.getEntries())
46 for mc in mcParticles:
47 p = mc.getMomentum()
48 h_momentum.Fill(p.Mag())
49 h_pt.Fill(p.Perp())
50 h_phi.Fill(p.Phi() / math.pi * 180)
51 h_theta.Fill(p.Theta() / math.pi * 180)
52 h_costheta.Fill(math.cos(p.Theta()))
53 h_pdg.Fill(mc.getPDG())
54 h_xyvertex.Fill(mc.getProductionVertex().X() * 1e4,
55 mc.getProductionVertex().Y() * 1e4)
56 h_zxvertex.Fill(mc.getProductionVertex().Z() * 1e4,
57 mc.getProductionVertex().X() * 1e4)
58 h_zyvertex.Fill(mc.getProductionVertex().Z() * 1e4,
59 mc.getProductionVertex().Y() * 1e4)
60
61
62# Create Modules
63eventinfosetter = register_module('EventInfoSetter')
64progress = register_module('Progress')
65particlegun = register_module('ParticleGun')
66vertexsmear = register_module('SmearPrimaryVertex')
67showMCPart = ShowMCParticles()
68
69# Specify number of events to generate
70eventinfosetter.param({'evtNumList': [10000], 'runList': [1]})
71
72# Set parameters for particlegun
73particlegun.param({
74 # Generate 5 tracks
75 'nTracks': 5,
76 # But vary the number of tracks according to Poisson distribution
77 'varyNTracks': True,
78 # Generate pi+, pi-, e+ and e-
79 'pdgCodes': [211, -211, 11, -11],
80 # having a normal distributed transversal momentum
81 'momentumGeneration': 'normalPt',
82 # with a mean of 5 GeV and a width of 1 GeV
83 'momentumParams': [5.0, 1.0],
84 # use normal distribution for phi
85 'phiGeneration': 'normal',
86 # with a mean of 90 degree and w width of 30 degree
87 'phiParams': [90, 30],
88 # and generate theta uniform in cos(theta)
89 'thetaGeneration': 'uniformCos',
90 # between 17 and 150 degree
91 'thetaParams': [17, 150],
92 # finally, vertex generation, fixed, we use smearing module
93 'vertexGeneration': 'fixed',
94 'xVertexParams': [0],
95 'yVertexParams': [0],
96 'zVertexParams': [0],
97 # all tracks should originate from the same vertex
98 'independentVertices': False,
99})
100
101
102# The beamspot size is determined by using the beamsize for LER and HER from
103# the TDR, multiply both zx-distributions including their crossing angle and
104# fit a rotated 2D gaussian to the resulting shape.
105# Parameters used (angle is the angle between beam and z axis):
106# HER: sigma_x = 10.2µm, sigma_z = 6mm, angle= 41.5mrad
107# LER: sigma_x = 7.76µm, sigma_z = 5mm, angle=-41.5mrad
108vertexsmear.param({
109 # Beamspot size in x
110 'sigma_pvx': 6.18e-4, # µm
111 # Beamspot size in y
112 'sigma_pvy': 59e-7, # nm
113 # Beamspot size in z
114 'sigma_pvz': 154e-4, # µm
115 # Angle between beamspot and z axis, rotation around y
116 'angle_pv_zx': -1.11e-2,
117})
118
119
120# create processing path
121main = create_path()
122main.add_module(eventinfosetter)
123main.add_module(progress)
124main.add_module(particlegun)
125main.add_module(vertexsmear)
126main.add_module(showMCPart)
127
128# generate events
129process(main)
130
131# show call statistics
132print(statistics)
133
134# Create a Canvas to show histograms
135c = TCanvas('Canvas', 'Canvas', 1920, 768)
136c.Divide(5, 2, 1e-5, 1e-5)
137
138# Draw all histograms
139histograms = [h_nTracks, h_pdg, h_momentum, h_pt, h_theta, h_costheta, h_phi]
140vertexhists = [h_xyvertex, h_zxvertex, h_zyvertex]
141for (i, h) in enumerate(histograms):
142 c.cd(i + 1)
143 h.SetMinimum(0)
144 h.Draw()
145for (j, h) in enumerate(vertexhists, i + 2):
146 c.cd(j)
147 h.Draw('colz')
148c.Update()
149
150# Wait for enter to be pressed
151print("Press Enter to exit ...")
152sys.stdin.readline()
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72