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