Belle II Software development
ClusterEfficiency Class Reference
Inheritance diagram for ClusterEfficiency:

Public Member Functions

def initialize (self)
 
def terminate (self)
 
def fill_truehits (self, phi, p, truehits)
 
def event (self)
 

Public Attributes

 rfile
 Output file to store all plots.
 
 profiles
 flat list of all profiles for easy access
 
 eff
 layer/momentum hierarchy of all profiles
 

Detailed Description

Plot Efficiency to find a U and V cluster for each given truehit.

Definition at line 36 of file ClusterEfficiency.py.

Member Function Documentation

◆ event()

def event (   self)
Update the efficiencies by iterating over all primary particles

Definition at line 151 of file ClusterEfficiency.py.

151 def event(self):
152 """
153 Update the efficiencies by iterating over all primary particles
154 """
155 mcparticles = Belle2.PyStoreArray("MCParticles")
156 for mcp in mcparticles:
157 # we only look at primary particles
158 if not mcp.hasStatus(1):
159 continue
160
161 # let's get the momentum and determine which of the ones we wanted
162 # to generate it is
163 p = mcp.getMomentum()
164 p_gen = None
165 for i in momenta:
166 if abs(p.R() - i) < 0.05:
167 p_gen = i
168 break
169
170 # meh, something strange with the momentum, ignore this one
171 if p_gen is None:
172 b2.B2WARNING(f"Strange particle momentum: {p.R():f}, expected one of {', '.join('' for p in momenta)}")
173 continue
174
175 # and check all truehits
176 pxdtruehits = mcp.getRelationsTo("PXDTrueHits")
177 self.fill_truehits(math.degrees(p.Phi()), p_gen, pxdtruehits)
178 svdtruehits = mcp.getRelationsTo("SVDTrueHits")
179 self.fill_truehits(math.degrees(p.Phi()), p_gen, svdtruehits)
180
181
182# Now let's create a path to simulate our events. We need a bit of statistics but
183# that's not too bad since we only simulate single muons
A (simplified) python wrapper for StoreArray.
Definition: PyStoreArray.h:72

◆ fill_truehits()

def fill_truehits (   self,
  phi,
  p,
  truehits 
)
Loop over all truehits for a track with the given angle phi and momentum
p and update the efficiency profiles.

Definition at line 110 of file ClusterEfficiency.py.

110 def fill_truehits(self, phi, p, truehits):
111 """
112 Loop over all truehits for a track with the given angle phi and momentum
113 p and update the efficiency profiles.
114 """
115 # iterate over all truehits and look for the clusters
116 for i, truehit in enumerate(truehits):
117 # we ignore all truehits which have an negative weight
118 if truehits.weight(i) < 0:
119 continue
120
121 # get layer number and a list of all clusters related to the truehit
122 layer = truehit.getSensorID().getLayerNumber()
123 if isinstance(truehit, Belle2.SVDTrueHit):
124 is_pxd = False
125 clusters = truehit.getRelationsFrom("SVDClusters")
126 else:
127 is_pxd = True
128 clusters = truehit.getRelationsFrom("PXDClusters")
129
130 # now check if we find a cluster
131 has_cluster = {False: 0, True: 0}
132 for j, cls in enumerate(clusters):
133 # we ignore all clusters where less then 100 electrons come from
134 # our truehit
135 if clusters.weight(j) < 100:
136 continue
137
138 if is_pxd:
139 # for pxt we always have both directions so set efficiency
140 # for this TrueHit to 1
141 has_cluster[True] = 1
142 has_cluster[False] = 1
143 else:
144 # for SVD we remember that we set the efficiency for the
145 # direction of the cluster
146 has_cluster[cls.isUCluster()] = 1
147
148 # and we fill the profile
149 self.eff[layer][p].Fill(phi, has_cluster[True] & has_cluster[False])
150
Class SVDTrueHit - Records of tracks that either enter or leave the sensitive volume.
Definition: SVDTrueHit.h:33

◆ initialize()

def initialize (   self)
Create ROOT TProfiles for all layers and momenta.

Definition at line 41 of file ClusterEfficiency.py.

41 def initialize(self):
42 """
43 Create ROOT TProfiles for all layers and momenta.
44 """
45
46 # let's create a root file to store all profiles
47
48 self.rfile = ROOT.TFile("ClusterEfficiency.root", "RECREATE")
49 self.rfile.cd()
50
51 self.profiles = []
52
53 self.eff = {}
54
55 def profile(name, title, text, contact):
56 """small helper function to create a phi profile and set the names
57 and descriptions"""
58 prof = ROOT.TProfile(name, title, 60, -180, 180)
59 prof.GetListOfFunctions().Add(ROOT.TNamed("Description", text))
60 prof.GetListOfFunctions().Add(ROOT.TNamed("Contact", contact))
61 prof.GetListOfFunctions().Add(ROOT.TNamed("Check", "Should be close to 1 everywhere"))
62 # make a list of all profiles
63 self.profiles.append(prof)
64 return prof
65
66 # name, title and description templates
67 prof_name = "ClusterEfficiency_layer{layer}_{p:.1f}GeV"
68 prof_title = "Cluster Efficiency in #phi, layer={layer}, "\
69 + "p={p:.1f} GeV;#phi in degrees;efficiency"
70 prof_text = "Efficiency to find a U and V cluster for any given truehit "\
71 + "in layer {layer} when simulation muons with p={p:.1f} GeV uniformly "\
72 + "in phi and theta={theta[0]}+-{theta[1]} degree. phi is the angle "\
73 + "of the generated particle, not the truehit position."
74 # contact person for PXD (layer<3 = True) or SVD (layer<3 = False)
75 prof_contact = {
76 True: "Benjamin Schwenker <Benjamin.Schwenker@phys.uni-goettingen.de>",
77 False: "Andrzej Bozek <bozek@belle2.ifj.edu.pl>",
78 }
79
80 # create all profiles
81 for layer in range(1, 7):
82 self.eff[layer] = {}
83 for p in momenta:
84 name = prof_name.format(p=p, layer=layer)
85 title = prof_title.format(p=p, layer=layer)
86 text = prof_text.format(p=p, layer=layer, theta=theta_params)
87 self.eff[layer][p] = profile(name, title, text, prof_contact[layer < 3])
88

◆ terminate()

def terminate (   self)
Format all profiles and write the ROOT file.

Definition at line 89 of file ClusterEfficiency.py.

89 def terminate(self):
90 """
91 Format all profiles and write the ROOT file.
92 """
93 # determine min value in all profiles to set all profiles to a common
94 # range. We always want to show at least the range from 0.9
95 minVal = 0.9
96 for prof in self.profiles:
97 minBin = prof.GetMinimumBin()
98 # minimum value is bin content - error and a 5% margin
99 minVal = min(minVal, (prof.GetBinContent(minBin) - prof.GetBinError(minBin)) * 0.95)
100
101 # no let's set the range of all profiles to be equal
102 for prof in self.profiles:
103 prof.SetMinimum(max(0, minVal))
104 prof.SetMaximum(1.02)
105 prof.SetStats(False)
106
107 self.rfile.Write()
108 self.rfile.Close()
109

Member Data Documentation

◆ eff

eff

layer/momentum hierarchy of all profiles

Definition at line 53 of file ClusterEfficiency.py.

◆ profiles

profiles

flat list of all profiles for easy access

Definition at line 51 of file ClusterEfficiency.py.

◆ rfile

rfile

Output file to store all plots.

Definition at line 48 of file ClusterEfficiency.py.


The documentation for this class was generated from the following file: