Belle II Software  release-08-01-10
ClusterEfficiency.py
1 #!/usr/bin/env python3
2 # -*- coding: utf-8 -*-
3 
4 
11 
12 """
13 <header>
14  <contact>Benjamin.Schwenker@phys.uni-goettingen.de</contact>
15  <description>
16  Plot the efficiency to find a cluster for a given truehit for all layers of
17  the VXD.
18  </description>
19 </header>
20 """
21 import math
22 import basf2 as b2
23 import svd
24 import pxd
25 import ROOT
26 from ROOT import Belle2
27 b2.set_log_level(b2.LogLevel.ERROR)
28 # set some random seed
29 b2.set_random_seed(10346)
30 
31 # momenta to generate the plots for
32 momenta = [3.0]
33 # theta parameters
34 theta_params = [90, 0.1]
35 
36 
37 class ClusterEfficiency(b2.Module):
38  """
39  Plot Efficiency to find a U and V cluster for each given truehit.
40  """
41 
42  def initialize(self):
43  """
44  Create ROOT TProfiles for all layers and momenta.
45  """
46 
47  # let's create a root file to store all profiles
48 
49  self.rfilerfile = ROOT.TFile("ClusterEfficiency.root", "RECREATE")
50  self.rfilerfile.cd()
51 
52  self.profilesprofiles = []
53 
54  self.effeff = {}
55 
56  def profile(name, title, text, contact):
57  """small helper function to create a phi profile and set the names
58  and descriptions"""
59  prof = ROOT.TProfile(name, title, 60, -180, 180)
60  prof.GetListOfFunctions().Add(ROOT.TNamed("Description", text))
61  prof.GetListOfFunctions().Add(ROOT.TNamed("Contact", contact))
62  prof.GetListOfFunctions().Add(ROOT.TNamed("Check", "Should be close to 1 everywhere"))
63  # make a list of all profiles
64  self.profilesprofiles.append(prof)
65  return prof
66 
67  # name, title and description templates
68  prof_name = "ClusterEfficiency_layer{layer}_{p:.1f}GeV"
69  prof_title = "Cluster Efficiency in #phi, layer={layer}, "\
70  + "p={p:.1f} GeV;#phi in degrees;efficiency"
71  prof_text = "Efficiency to find a U and V cluster for any given truehit "\
72  + "in layer {layer} when simulation muons with p={p:.1f} GeV uniformly "\
73  + "in phi and theta={theta[0]}+-{theta[1]} degree. phi is the angle "\
74  + "of the generated particle, not the truehit position."
75  # contact person for PXD (layer<3 = True) or SVD (layer<3 = False)
76  prof_contact = {
77  True: "Benjamin Schwenker <Benjamin.Schwenker@phys.uni-goettingen.de>",
78  False: "Andrzej Bozek <bozek@belle2.ifj.edu.pl>",
79  }
80 
81  # create all profiles
82  for layer in range(1, 7):
83  self.effeff[layer] = {}
84  for p in momenta:
85  name = prof_name.format(p=p, layer=layer)
86  title = prof_title.format(p=p, layer=layer)
87  text = prof_text.format(p=p, layer=layer, theta=theta_params)
88  self.effeff[layer][p] = profile(name, title, text, prof_contact[layer < 3])
89 
90  def terminate(self):
91  """
92  Format all profiles and write the ROOT file.
93  """
94  # determine min value in all profiles to set all profiles to a common
95  # range. We always want to show at least the range from 0.9
96  minVal = 0.9
97  for prof in self.profilesprofiles:
98  minBin = prof.GetMinimumBin()
99  # minimum value is bin content - error and a 5% margin
100  minVal = min(minVal, (prof.GetBinContent(minBin) - prof.GetBinError(minBin)) * 0.95)
101 
102  # no let's set the range of all profiles to be equal
103  for prof in self.profilesprofiles:
104  prof.SetMinimum(max(0, minVal))
105  prof.SetMaximum(1.02)
106  prof.SetStats(False)
107 
108  self.rfilerfile.Write()
109  self.rfilerfile.Close()
110 
111  def fill_truehits(self, phi, p, truehits):
112  """
113  Loop over all truehits for a track with the given angle phi and momentum
114  p and update the efficiency profiles.
115  """
116  # iterate over all truehits and look for the clusters
117  for i, truehit in enumerate(truehits):
118  # we ignore all truehits which have an negative weight
119  if truehits.weight(i) < 0:
120  continue
121 
122  # get layer number and a list of all clusters related to the truehit
123  layer = truehit.getSensorID().getLayerNumber()
124  if isinstance(truehit, Belle2.SVDTrueHit):
125  is_pxd = False
126  clusters = truehit.getRelationsFrom("SVDClusters")
127  else:
128  is_pxd = True
129  clusters = truehit.getRelationsFrom("PXDClusters")
130 
131  # now check if we find a cluster
132  has_cluster = {False: 0, True: 0}
133  for j, cls in enumerate(clusters):
134  # we ignore all clusters where less then 100 electrons come from
135  # our truehit
136  if clusters.weight(j) < 100:
137  continue
138 
139  if is_pxd:
140  # for pxt we always have both directions so set efficiency
141  # for this TrueHit to 1
142  has_cluster[True] = 1
143  has_cluster[False] = 1
144  else:
145  # for SVD we remember that we set the efficiency for the
146  # direction of the cluster
147  has_cluster[cls.isUCluster()] = 1
148 
149  # and we fill the profile
150  self.effeff[layer][p].Fill(phi, has_cluster[True] & has_cluster[False])
151 
152  def event(self):
153  """
154  Update the efficiencies by iterating over all primary particles
155  """
156  mcparticles = Belle2.PyStoreArray("MCParticles")
157  for mcp in mcparticles:
158  # we only look at primary particles
159  if not mcp.hasStatus(1):
160  continue
161 
162  # let's get the momentum and determine which of the ones we wanted
163  # to generate it is
164  p = mcp.getMomentum()
165  p_gen = None
166  for i in momenta:
167  if abs(p.R() - i) < 0.05:
168  p_gen = i
169  break
170 
171  # meh, something strange with the momentum, ignore this one
172  if p_gen is None:
173  b2.B2WARNING("Strange particle momentum: %f, expected one of %s" %
174  (p.R(), ", ".join(str() for p in momenta)))
175  continue
176 
177  # and check all truehits
178  pxdtruehits = mcp.getRelationsTo("PXDTrueHits")
179  self.fill_truehitsfill_truehits(math.degrees(p.Phi()), p_gen, pxdtruehits)
180  svdtruehits = mcp.getRelationsTo("SVDTrueHits")
181  self.fill_truehitsfill_truehits(math.degrees(p.Phi()), p_gen, svdtruehits)
182 
183 
184 # Now let's create a path to simulate our events. We need a bit of statistics but
185 # that's not too bad since we only simulate single muons
186 main = b2.create_path()
187 main.add_module("EventInfoSetter", evtNumList=[10000])
188 main.add_module("Gearbox")
189 # we only need the vxd for this
190 main.add_module("Geometry", components=['MagneticFieldConstant4LimitedRSVD',
191  'BeamPipe', 'PXD', 'SVD'])
192 particlegun = main.add_module("ParticleGun")
193 particlegun.param({
194  "nTracks": 1,
195  "pdgCodes": [13, -13],
196  # generate discrete momenta with equal weights
197  "momentumGeneration": 'discrete',
198  "momentumParams": momenta + [1]*len(momenta),
199  "thetaGeneration": 'normal',
200  "thetaParams": theta_params,
201 })
202 main.add_module("FullSim")
207 
208 clusterefficiency = ClusterEfficiency()
209 main.add_module(clusterefficiency)
210 main.add_module("Progress")
211 
212 b2.process(main)
213 print(b2.statistics)
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72
Class SVDTrueHit - Records of tracks that either enter or leave the sensitive volume.
Definition: SVDTrueHit.h:33
eff
layer/momentum hierarchy of all profiles
rfile
Output file to store all plots.
profiles
flat list of all profiles for easy access
def fill_truehits(self, phi, p, truehits)
def add_pxd_reconstruction(path, clusterName=None, digitsName=None, usePXDClusterShapes=False, spacePointsName='PXDSpacePoints')
Definition: __init__.py:105
def add_pxd_simulation(path, digitsName=None, activatePixelMasks=True, activateGainCorrection=True)
Definition: __init__.py:147
def add_svd_simulation(path, useConfigFromDB=False, daqMode=2, relativeShift=9)
Definition: __init__.py:245
def add_svd_reconstruction(path, isROIsimulation=False, createRecoDigits=False, applyMasking=False)
Definition: __init__.py:15