Belle II Software  release-08-01-10
extmuid_createPlots.py
1 #!/usr/bin/env python3
2 
3 
10 
11 
17 
18 """
19 <header>
20  <input>muon-ExtMuidValidation.root,
21  pion-ExtMuidValidation.root</input>
22  <contact>giacomo.pietro@kit.edu</contact>
23  <description>Create validation plots for Ext and Muid</description>
24 </header>
25 """
26 
27 import ROOT
28 
29 ROOT.PyConfig.IgnoreCommandLineOptions = True # noqa
30 ROOT.gROOT.SetBatch(True) # noqa
31 
32 from ROOT import TChain, TFile, TH1F, TH2F, TNamed, gStyle
33 import sys
34 import math
35 import numpy as np
36 from optparse import OptionParser
37 
38 # contact person information
39 # is added to the plot descriptions
40 CONTACT_PERSON = {'Name': 'Giacomo De Pietro',
41  'Email': 'giacomo.pietro@kit.edu'}
42 
43 ACTIVE = True
44 
45 
46 def main():
47  """Create validation plots for Ext and Muid"""
48 
49  option_parser = OptionParser()
50  option_parser.add_option('-i', '--input-file', dest='input_file',
51  default='../*on-ExtMuidValidation.root',
52  help='Root file with Ext and Muid validation data.'
53  )
54  option_parser.add_option('-o', '--output-file', dest='output_file',
55  default='ExtMuidValidation.root',
56  help='Root file with Ext and Muid validation histograms.')
57 
58  (options, args) = option_parser.parse_args()
59 
60  # load chain of input files
61  file_chain = TChain('tree')
62  file_chain.Add(options.input_file)
63 
64  # for file in file_chain:
65 
66  number_entries = 0
67  try:
68  number_entries = file_chain.GetEntries()
69  except AttributeError:
70  print('Could not load input file(s) %s.' % options.input_file)
71 
72  if number_entries == 0:
73  print('Data tree is empty or does not exist in file(s) %s. Exit.' % options.input_file)
74  sys.exit(0)
75 
76  # open the output root file
77  output_root_file = TFile(options.output_file, 'recreate')
78 
79  # create and draw histograms
80  gStyle.SetOptStat(0)
81  draw_exthits(file_chain)
82  draw_likelihoods(file_chain)
83 
84  # close the output file
85  output_root_file.Write()
86  output_root_file.Close()
87 
88 
89 def draw_exthits(file_chain):
90  """
91  Draw the ExtHit-related distributions.
92  """
93 
94  contact_mail = "giacomo.pietro@kit.edu"
95 
96  # NOTE: *.Draw() must precede *.GetListOfFunctions().Add() or the latter will be discarded!
97  detectorID = TH1F('DetectorID', 'Detector ID for ExtHits', 8, -0.5, 7.5)
98  file_chain.Draw('ExtHits.m_DetectorID&0x0F>>DetectorID', '')
99  detectorID.GetXaxis().SetTitle('0=undefined, 1=PXD, 2=SVD, 3=CDC, 4=TOP, 5=ARICH, 6=ECL, 7=KLM')
100  detectorID.GetListOfFunctions().Add(TNamed('Description', "0=undefined, 1=PXD, 2=SVD, 3=CDC, 4=TOP, 5=ARICH, 6=ECL, 7=KLM"))
101  detectorID.GetListOfFunctions().Add(TNamed('Check', "ECL > (KLM ~ TOP) >> ARICH; others ~ 0"))
102  detectorID.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
103  detectorID.GetListOfFunctions().Add(TNamed('MetaOptions', 'shifter,pvalue-warn=0.50,pvalue-error=0.10'))
104  detectorID.Write()
105 
106  tof = TH1F('TOF', 'Time of Flight for non-KLM ExtHits', 100, 0.0, 25.0)
107  file_chain.Draw('ExtHits.m_TOF>>TOF', '(ExtHits.m_DetectorID&0x0F)!=7')
108  tof.GetXaxis().SetTitle('t (ns)')
109  tof.GetListOfFunctions().Add(TNamed('Description', "Time of flight along track from IP to each ExtHit"))
110  tof.GetListOfFunctions().Add(TNamed('Check', "Dominant peak at ~4.5 ns, secondary peak at ~6 ns"))
111  tof.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
112  tof.GetListOfFunctions().Add(TNamed('MetaOptions', 'expert,pvalue-warn=0.50,pvalue-error=0.10'))
113  tof.Write()
114 
115  tofKLM = TH1F('TOFKLM', 'Time of Flight for KLM ExtHits', 100, 0.0, 25.0)
116  file_chain.Draw('ExtHits.m_TOF>>TOFKLM', '(ExtHits.m_DetectorID&0x0F)==7')
117  tofKLM.GetXaxis().SetTitle('t (ns)')
118  tofKLM.GetListOfFunctions().Add(TNamed('Description', "Time of flight along track from IP to each ExtHit"))
119  tofKLM.GetListOfFunctions().Add(TNamed('Check', "Broad rising distribution at 7-11 ns, then falling to ~17 ns"))
120  tofKLM.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
121  tofKLM.GetListOfFunctions().Add(TNamed('MetaOptions', 'expert,pvalue-warn=0.50,pvalue-error=0.10'))
122  tofKLM.Write()
123 
124  r = TH1F('r', 'r for non-KLM ExtHits', 40, 0.0, 200.0)
125  file_chain.Draw('ExtHits.getPosition().Perp()>>r', '(ExtHits.m_DetectorID&0x0F)!=7')
126  r.GetXaxis().SetTitle('r (cm)')
127  r.GetListOfFunctions().Add(TNamed('Description', "Radial position in transverse plane of each ExtHit"))
128  r.GetListOfFunctions().Add(TNamed('Check', "Sharp peak at 120 cm; broad peak at 140 cm"))
129  r.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
130  r.GetListOfFunctions().Add(TNamed('MetaOptions', 'expert,pvalue-warn=0.50,pvalue-error=0.10'))
131  r.Write()
132 
133  z = TH1F('z', 'z for non-KLM ExtHits', 100, -200.0, 300.0)
134  file_chain.Draw('ExtHits.getPosition().Z()>>z', '(ExtHits.m_DetectorID&0x0F)!=7')
135  z.GetXaxis().SetTitle('z (cm)')
136  z.GetListOfFunctions().Add(TNamed('Description', "Axial position of each ExtHit"))
137  z.GetListOfFunctions().Add(TNamed('Check', "Broad peak centered at 0 cm; sharp dip at 0; secondary peak at 170"))
138  z.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
139  z.GetListOfFunctions().Add(TNamed('MetaOptions', 'expert,pvalue-warn=0.50,pvalue-error=0.10'))
140  z.Write()
141 
142  rKLM = TH1F('rKLM', 'r for KLM ExtHits', 50, 100.0, 350.0)
143  file_chain.Draw('ExtHits.getPosition().Perp()>>rKLM', '(ExtHits.m_DetectorID&0x0F)==7')
144  rKLM.GetXaxis().SetTitle('r (cm)')
145  rKLM.GetListOfFunctions().Add(TNamed('Description', "Radial position in transverse plane of each ExtHit"))
146  rKLM.GetListOfFunctions().Add(TNamed('Check', "Low shoulder below 200 cm (EKLM); comb-like pattern above 200 cm (BKLM)"))
147  rKLM.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
148  rKLM.GetListOfFunctions().Add(TNamed('MetaOptions', 'expert,pvalue-warn=0.50,pvalue-error=0.10'))
149  rKLM.Write()
150 
151  zKLM = TH1F('zKLM', 'z for KLM ExtHits', 140, -300.0, 400.0)
152  file_chain.Draw('ExtHits.getPosition().Z()>>zKLM', '(ExtHits.m_DetectorID&0x0F)==7')
153  zKLM.GetXaxis().SetTitle('z (cm)')
154  zKLM.GetListOfFunctions().Add(TNamed('Description', "Axial position of each ExtHit"))
155  zKLM.GetListOfFunctions().Add(TNamed('Check', "Comb-like pattern at z<-180 and z>280 cm (EKLM); broad peak at 0 (BKLM)"))
156  zKLM.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
157  zKLM.GetListOfFunctions().Add(TNamed('MetaOptions', 'expert,pvalue-warn=0.50,pvalue-error=0.10'))
158  zKLM.Write()
159 
160  xy = TH2F('xy', 'y vs x for non-KLM ExtHits', 100, -200.0, 200.0, 100, -200.0, 200.0)
161  file_chain.Draw('ExtHits.getPosition().Y():ExtHits.getPosition().X()>>xy', '(ExtHits.m_DetectorID&0x0F)!=7')
162  xy.GetXaxis().SetTitle('x (cm)')
163  xy.GetYaxis().SetTitle('y (cm)')
164  xy.GetListOfFunctions().Add(TNamed('Description', "Position projected into transverse plane of each ExtHit"))
165  xy.GetListOfFunctions().Add(TNamed('Check', "Cylindrical, with most hits in ECL and TOP"))
166  xy.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
167  xy.GetListOfFunctions().Add(TNamed('MetaOptions', "box, expert"))
168  xy.Write()
169 
170  xz = TH2F('xz', 'x vs z for non-KLM ExtHits', 125, -200.0, 300.0, 100, -200.0, 200.0)
171  file_chain.Draw('ExtHits.getPosition().X():ExtHits.getPosition().Z()>>xz', '(ExtHits.m_DetectorID&0x0F)!=7')
172  xz.GetXaxis().SetTitle('z (cm)')
173  xz.GetYaxis().SetTitle('x (cm)')
174  xz.GetListOfFunctions().Add(TNamed('Description', "Position projected into x-z plane of each ExtHit"))
175  xz.GetListOfFunctions().Add(TNamed('Check', " "))
176  xz.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
177  xz.GetListOfFunctions().Add(TNamed('MetaOptions', "box, expert"))
178  xz.Write()
179 
180  yz = TH2F('yz', 'y vs z for non-KLM ExtHits', 125, -200.0, 300.0, 100, -200.0, 200.0)
181  file_chain.Draw('ExtHits.getPosition().Y():ExtHits.getPosition().Z()>>yz', '(ExtHits.m_DetectorID&0x0F)!=7')
182  yz.GetXaxis().SetTitle('z (cm)')
183  yz.GetYaxis().SetTitle('y (cm)')
184  yz.GetListOfFunctions().Add(TNamed('Description', "Position projected into y-z plane of each ExtHit"))
185  yz.GetListOfFunctions().Add(TNamed('Check', " "))
186  yz.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
187  yz.GetListOfFunctions().Add(TNamed('MetaOptions', "box, expert"))
188  yz.Write()
189 
190  xyKLM = TH2F('xyKLM', 'y vs x for KLM ExtHits', 140, -350.0, 350.0, 140, -350.0, 350.0)
191  file_chain.Draw('ExtHits.getPosition().Y():ExtHits.getPosition().X()>>xyKLM', '(ExtHits.m_DetectorID&0x0F)==7')
192  xyKLM.GetXaxis().SetTitle('x (cm)')
193  xyKLM.GetYaxis().SetTitle('y (cm)')
194  xyKLM.GetListOfFunctions().Add(TNamed('Description', "Position projected into transverse plane of each ExtHit"))
195  xyKLM.GetListOfFunctions().Add(TNamed('Check', "Octagonal (BKLM) and cylindrical (EKLM)"))
196  xyKLM.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
197  xyKLM.GetListOfFunctions().Add(TNamed('MetaOptions', "box, expert"))
198  xyKLM.Write()
199 
200  xzKLM = TH2F('xzKLM', 'x vs z for KLM ExtHits', 140, -300.0, 400.0, 140, -350.0, 350.0)
201  file_chain.Draw('ExtHits.getPosition().X():ExtHits.getPosition().Z()>>xzKLM', '(ExtHits.m_DetectorID&0x0F)==7')
202  xzKLM.GetXaxis().SetTitle('z (cm)')
203  xzKLM.GetYaxis().SetTitle('x (cm)')
204  xzKLM.GetListOfFunctions().Add(TNamed('Description', "Position projected into x-z plane of each ExtHit"))
205  xzKLM.GetListOfFunctions().Add(TNamed('Check', " "))
206  xzKLM.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
207  xzKLM.GetListOfFunctions().Add(TNamed('MetaOptions', "box, expert"))
208  xzKLM.Write()
209 
210  yzKLM = TH2F('yzKLM', 'y vs z for KLM ExtHits', 140, -300.0, 400.0, 140, -350.0, 350.0)
211  file_chain.Draw('ExtHits.getPosition().Y():ExtHits.getPosition().Z()>>yzKLM', '(ExtHits.m_DetectorID&0x0F)==7')
212  yzKLM.GetXaxis().SetTitle('z (cm)')
213  yzKLM.GetYaxis().SetTitle('y (cm)')
214  yzKLM.GetListOfFunctions().Add(TNamed('Description', "Position projected into y-z plane of each ExtHit"))
215  yzKLM.GetListOfFunctions().Add(TNamed('Check', " "))
216  yzKLM.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
217  yzKLM.GetListOfFunctions().Add(TNamed('MetaOptions', "box, expert"))
218  yzKLM.Write()
219 
220 
221 def draw_likelihoods(file_chain):
222  """
223  Draw the Muid likelihood-based distributions.
224  """
225 
226  contact_mail = "giacomo.pietro@kit.edu"
227 
228  outcome = TH1F('Outcome', 'Outcome', 5, -0.5, 4.5)
229 
230  llBins = 120
231  llMin = -60.0
232  llMax = 0.0
233  llMu_mu = TH1F('LLmu-muons', 'logL(mu) for true muons', llBins, llMin, llMax)
234  llMu_pi = TH1F('LLmu-pions', 'logL(mu) for true pions', llBins, llMin, llMax)
235  llPi_mu = TH1F('LLpi-muons', 'logL(pi) for true muons', llBins, llMin, llMax)
236  llPi_pi = TH1F('LLpi-pions', 'logL(pi) for true pions', llBins, llMin, llMax)
237  llDiff_mu = TH1F('LLdiff-muons', 'logL(mu)-logL(pi) for true muons', llBins * 2, llMin, -llMin)
238  llDiff_pi = TH1F('LLdiff-pions', 'logL(mu)-logL(pi) for true pions', llBins * 2, llMin, -llMin)
239 
240  layerBins = 15
241  layerMin = -0.5
242  layerMax = 14.5
243  barrelLayer_mu = TH1F('BarrelLayer-muons', 'Last barrel layer with hit for true muons', layerBins, layerMin, layerMax)
244  barrelLayer_pi = TH1F('BarrelLayer-pions', 'Last barrel layer with hit for true pions', layerBins, layerMin, layerMax)
245  endcapLayer_mu = TH1F('EndcapLayer-muons', 'Last endcap layer with hit for true muons', layerBins, layerMin, layerMax)
246  endcapLayer_pi = TH1F('EndcapLayer-pions', 'Last endcap layer with hit for true pions', layerBins, layerMin, layerMax)
247 
248  diffLayerBins = 20
249  diffLayerMin = -0.5
250  diffLayerMax = 19.5
251  diffLayer_mu = TH1F('DiffLayer-muons', 'Range difference for true muons', diffLayerBins, diffLayerMin, diffLayerMax)
252  diffLayer_pi = TH1F('DiffLayer-pions', 'Range difference for true pions', diffLayerBins, diffLayerMin, diffLayerMax)
253 
254  ndofBins = 32
255  ndofMin = -0.5
256  ndofMax = 31.5
257  ndof_mu = TH1F('Ndof-muons', 'Degrees of freedom for true muons', ndofBins, ndofMin, ndofMax)
258  ndof_pi = TH1F('Ndof-pions', 'Degrees of freedom for true pions', ndofBins, ndofMin, ndofMax)
259 
260  rchisqBins = 50
261  rchisqMin = 0.0
262  rchisqMax = 10.0
263  rchisq_mu = TH1F('rchisq-muons', 'Reduced chi-squared for true muons', rchisqBins, rchisqMin, rchisqMax)
264  rchisq_pi = TH1F('rchisq-pions', 'Reduced chi-squared for true pions', rchisqBins, rchisqMin, rchisqMax)
265 
266  momentumBins = 25 # for efficiency
267  momentumBins2 = 5 # for fake rate (low statistics)
268  momentumMin = 0.0
269  momentumMax = 5.0
270  efficiency_momentum = TH1F('Eff-momentum', 'Muon efficiency vs momentum', momentumBins, momentumMin, momentumMax)
271  efficiency_momentum_denom = TH1F('Eff-momentum-denom', 'Muon efficiency vs momentum', momentumBins, momentumMin, momentumMax)
272  fakerate_momentum = TH1F('FakeRate-momentum', 'Pion fake rate vs momentum', momentumBins2, momentumMin, momentumMax)
273  fakerate_momentum_denom = TH1F('FakeRate-momentum-denom', 'Pion fake rate vs momentum', momentumBins2, momentumMin, momentumMax)
274 
275  thetaBins = 35 # for efficiency
276  thetaBins2 = 7 # for fake rate (low statistics)
277  thetaMin = 10.0
278  thetaMax = 150.0
279  efficiency_theta = TH1F('Eff-theta', 'Muon efficiency vs theta', thetaBins, thetaMin, thetaMax)
280  efficiency_theta_denom = TH1F('Eff-theta-denom', 'Muon efficiency vs theta', thetaBins, thetaMin, thetaMax)
281  fakerate_theta = TH1F('FakeRate-theta', 'Pion fake rate vs theta', thetaBins2, thetaMin, thetaMax)
282  fakerate_theta_denom = TH1F('FakeRate-theta-denom', 'Pion fake rate vs theta', thetaBins2, thetaMin, thetaMax)
283 
284  phiBins = 36 # for efficiency
285  phiBins2 = 8 # for fake rate (low statistics)
286  phiMin = 0.0
287  phiMax = 360.0
288  efficiency_phi = TH1F('Eff-phi', 'Muon efficiency vs phi', phiBins, phiMin, phiMax)
289  efficiency_phi_denom = TH1F('Eff-phi-denom', 'Muon efficiency vs phi', phiBins, phiMin, phiMax)
290  fakerate_phi = TH1F('FakeRate-phi', 'Pion fake rate vs phi', phiBins2, phiMin, phiMax)
291  fakerate_phi_denom = TH1F('FakeRate-phi-denom', 'Pion fake rate vs phi', phiBins2, phiMin, phiMax)
292 
293  for entry in file_chain:
294  mcps = entry.MCParticles
295  momentum = mcps[0].getMomentum()
296  muids = entry.KLMMuidLikelihoods
297  for i in range(muids.GetEntriesFast()):
298  outcome.Fill(muids[i].getOutcome())
299  if muids[i].getOutcome() > 0:
300  llMu = muids[i].getLogL_mu()
301  llPi = muids[i].getLogL_pi()
302  blayer = muids[i].getBarrelHitLayer()
303  elayer = muids[i].getEndcapHitLayer()
304  diffLayer = muids[i].getExtLayer() - muids[i].getHitLayer()
305  ndof = muids[i].getDegreesOfFreedom()
306  chisq = muids[i].getChiSquared()
307  rchisq = -1.0
308  if ndof > 0:
309  rchisq = chisq / ndof
310  p = momentum.R()
311  theta = momentum.Theta() * 180.0 / np.pi
312  phi = momentum.Phi() * 180.0 / np.pi
313  if phi < 0.0:
314  phi = phi + 360.0
315  if 'muon' in file_chain.GetCurrentFile().GetName():
316  llMu_mu.Fill(llMu)
317  llPi_mu.Fill(llPi)
318  llDiff_mu.Fill(llMu - llPi)
319  barrelLayer_mu.Fill(blayer)
320  endcapLayer_mu.Fill(elayer)
321  diffLayer_mu.Fill(diffLayer)
322  ndof_mu.Fill(ndof)
323  rchisq_mu.Fill(rchisq)
324  efficiency_momentum_denom.Fill(p)
325  efficiency_theta_denom.Fill(theta)
326  efficiency_phi_denom.Fill(phi)
327  if llMu > llPi:
328  efficiency_momentum.Fill(p)
329  efficiency_theta.Fill(theta)
330  efficiency_phi.Fill(phi)
331  else:
332  llMu_pi.Fill(llMu)
333  llPi_pi.Fill(llPi)
334  llDiff_pi.Fill(llMu - llPi)
335  barrelLayer_pi.Fill(blayer)
336  endcapLayer_pi.Fill(elayer)
337  diffLayer_pi.Fill(diffLayer)
338  ndof_pi.Fill(ndof)
339  rchisq_pi.Fill(rchisq)
340  fakerate_momentum_denom.Fill(p)
341  fakerate_theta_denom.Fill(theta)
342  fakerate_phi_denom.Fill(phi)
343  if llMu > llPi:
344  fakerate_momentum.Fill(p)
345  fakerate_theta.Fill(theta)
346  fakerate_phi.Fill(phi)
347 
348  for j in range(efficiency_momentum_denom.GetNbinsX()):
349  num = efficiency_momentum.GetBinContent(j + 1)
350  denom = efficiency_momentum_denom.GetBinContent(j + 1)
351  if denom > 0:
352  efficiency_momentum.SetBinContent(j + 1, num / denom)
353  efficiency_momentum.SetBinError(j + 1, math.sqrt(num * (denom - num) / (denom * denom * denom)))
354  for j in range(efficiency_theta_denom.GetNbinsX()):
355  num = efficiency_theta.GetBinContent(j + 1)
356  denom = efficiency_theta_denom.GetBinContent(j + 1)
357  if denom > 0:
358  efficiency_theta.SetBinContent(j + 1, num / denom)
359  efficiency_theta.SetBinError(j + 1, math.sqrt(num * (denom - num) / (denom * denom * denom)))
360  for j in range(efficiency_phi_denom.GetNbinsX()):
361  num = efficiency_phi.GetBinContent(j + 1)
362  denom = efficiency_phi_denom.GetBinContent(j + 1)
363  if denom > 0:
364  efficiency_phi.SetBinContent(j + 1, num / denom)
365  efficiency_phi.SetBinError(j + 1, math.sqrt(num * (denom - num) / (denom * denom * denom)))
366  for j in range(fakerate_momentum_denom.GetNbinsX()):
367  num = fakerate_momentum.GetBinContent(j + 1)
368  denom = fakerate_momentum_denom.GetBinContent(j + 1)
369  if denom > 0:
370  fakerate_momentum.SetBinContent(j + 1, num / denom)
371  fakerate_momentum.SetBinError(j + 1, math.sqrt(num * (denom - num) / (denom * denom * denom)))
372  for j in range(fakerate_theta_denom.GetNbinsX()):
373  num = fakerate_theta.GetBinContent(j + 1)
374  denom = fakerate_theta_denom.GetBinContent(j + 1)
375  if denom > 0:
376  fakerate_theta.SetBinContent(j + 1, num / denom)
377  fakerate_theta.SetBinError(j + 1, math.sqrt(num * (denom - num) / (denom * denom * denom)))
378  for j in range(fakerate_phi_denom.GetNbinsX()):
379  num = fakerate_phi.GetBinContent(j + 1)
380  denom = fakerate_phi_denom.GetBinContent(j + 1)
381  if denom > 0:
382  fakerate_phi.SetBinContent(j + 1, num / denom)
383  fakerate_phi.SetBinError(j + 1, math.sqrt(num * (denom - num) / (denom * denom * denom)))
384 
385  # NOTE: *.Fill() must precede *.GetListOfFunctions().Add() or the latter will be discarded!
386  outcome.GetListOfFunctions().Add(TNamed('Description', "0=not in KLM, 1/2=barrel/endcap stop, 3/4=barrel/endcap exit"))
387  outcome.GetListOfFunctions().Add(TNamed('Check', "Peak at 3; rest ~ flat"))
388  outcome.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
389  outcome.GetListOfFunctions().Add(TNamed('MetaOptions', 'expert,pvalue-warn=0.50,pvalue-error=0.10'))
390  outcome.SetMinimum(0.0)
391 
392  llMu_mu.GetListOfFunctions().Add(TNamed('Description', "Log-likelihood of muon hypothesis for true muons"))
393  llMu_mu.GetListOfFunctions().Add(TNamed('Check', "Sharp peak at 0 with tail on negative side"))
394  llMu_mu.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
395  llMu_mu.GetListOfFunctions().Add(TNamed('MetaOptions', 'expert,pvalue-warn=0.50,pvalue-error=0.10'))
396  llMu_pi.GetListOfFunctions().Add(TNamed('Description', "Log-likelihood of muon hypothesis for true pions"))
397  llMu_pi.GetListOfFunctions().Add(TNamed('Check', "Sharp double peak at ~52 with tail up to 0"))
398  llMu_pi.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
399  llMu_pi.GetListOfFunctions().Add(TNamed('MetaOptions', 'expert,pvalue-warn=0.50,pvalue-error=0.10'))
400  llPi_mu.GetListOfFunctions().Add(TNamed('Description', "Log-likelihood of pion hypothesis for true muons"))
401  llPi_mu.GetListOfFunctions().Add(TNamed('Check', "Sharp peak at ~-42 with tail up to 0"))
402  llPi_mu.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
403  llPi_mu.GetListOfFunctions().Add(TNamed('MetaOptions', 'expert,pvalue-warn=0.50,pvalue-error=0.10'))
404  llPi_pi.GetListOfFunctions().Add(TNamed('Description', "Log-likelihood of pion hypothesis for true pions"))
405  llPi_pi.GetListOfFunctions().Add(TNamed('Check', "Sharp peak at 0 with tail on negative side"))
406  llPi_pi.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
407  llPi_pi.GetListOfFunctions().Add(TNamed('MetaOptions', 'expert,pvalue-warn=0.50,pvalue-error=0.10'))
408  llDiff_mu.GetListOfFunctions().Add(TNamed('Description', "Log-likelihood difference of muon-pion hypotheses for true muons"))
409  llDiff_mu.GetListOfFunctions().Add(TNamed('Check', "Sharp peak at ~+42 with tail on lower side"))
410  llDiff_mu.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
411  llDiff_mu.GetListOfFunctions().Add(TNamed('MetaOptions', 'shifter,pvalue-warn=0.50,pvalue-error=0.10'))
412  llDiff_pi.GetListOfFunctions().Add(TNamed('Description', "Log-likelihood difference of muon-pion hypotheses for true pions"))
413  llDiff_pi.GetListOfFunctions().Add(TNamed('Check', "Sharp peak at ~-52 with tail on higher side"))
414  llDiff_pi.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
415  llDiff_pi.GetListOfFunctions().Add(TNamed('MetaOptions', 'shifter,pvalue-warn=0.50,pvalue-error=0.10'))
416 
417  barrelLayer_mu.GetListOfFunctions().Add(TNamed('Description', "Outermost barrel layer with matching hit for true muons"))
418  barrelLayer_mu.GetListOfFunctions().Add(TNamed('Check', "Sharp peak at layer 14"))
419  barrelLayer_mu.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
420  barrelLayer_mu.GetListOfFunctions().Add(TNamed('MetaOptions', 'expert,pvalue-warn=0.50,pvalue-error=0.10'))
421  barrelLayer_pi.GetListOfFunctions().Add(TNamed('Description', "Outermost barrel layer with matching hit for true pions"))
422  barrelLayer_pi.GetListOfFunctions().Add(TNamed('Check', "Peak at layers 0-1 with tail above that"))
423  barrelLayer_pi.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
424  barrelLayer_pi.GetListOfFunctions().Add(TNamed('MetaOptions', 'expert,pvalue-warn=0.50,pvalue-error=0.10'))
425  endcapLayer_mu.GetListOfFunctions().Add(TNamed('Description', "Outermost endcap layer with matching hit for true muons"))
426  endcapLayer_mu.GetListOfFunctions().Add(TNamed('Check', "Sharp peaks at layers 11 (backward) and 13 (forward)"))
427  endcapLayer_mu.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
428  endcapLayer_mu.GetListOfFunctions().Add(TNamed('MetaOptions', 'expert,pvalue-warn=0.50,pvalue-error=0.10'))
429  endcapLayer_pi.GetListOfFunctions().Add(TNamed('Description', "Outermost endcap layer with matching hit for true pions"))
430  endcapLayer_pi.GetListOfFunctions().Add(TNamed('Check', "Peak at layers 0-1 with tail above that"))
431  endcapLayer_pi.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
432  endcapLayer_pi.GetListOfFunctions().Add(TNamed('MetaOptions', 'expert,pvalue-warn=0.50,pvalue-error=0.10'))
433 
434  diffLayer_mu.GetListOfFunctions().Add(TNamed('Description', "Outermost extrapolated-vs-hit layer difference for true muons"))
435  diffLayer_mu.GetListOfFunctions().Add(TNamed('Check', "Sharp peak at 0"))
436  diffLayer_mu.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
437  diffLayer_mu.GetListOfFunctions().Add(TNamed('MetaOptions', 'shifter,pvalue-warn=0.50,pvalue-error=0.10'))
438  diffLayer_mu.GetXaxis().SetTitle('Last Ext Layer - Last Hit Layer')
439  diffLayer_pi.GetListOfFunctions().Add(TNamed('Description', "Outermost extrapolated-vs-hit layer difference for true pions"))
440  diffLayer_pi.GetListOfFunctions().Add(TNamed('Check', "Sharp peak at 15 with tail below that"))
441  diffLayer_pi.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
442  diffLayer_pi.GetListOfFunctions().Add(TNamed('MetaOptions', 'shifter,pvalue-warn=0.50,pvalue-error=0.10'))
443  diffLayer_pi.GetXaxis().SetTitle('Last Ext Layer - Last Hit Layer')
444 
445  ndof_mu.GetListOfFunctions().Add(TNamed('Description', "(Number of matching-hit layers) x 2 for true muons"))
446  ndof_mu.GetListOfFunctions().Add(TNamed('Check', "Even-only entries; peak at 30 with tail below"))
447  ndof_mu.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
448  ndof_mu.GetListOfFunctions().Add(TNamed('MetaOptions', 'expert,pvalue-warn=0.50,pvalue-error=0.10'))
449  ndof_mu.GetXaxis().SetTitle('#dof')
450  ndof_pi.GetListOfFunctions().Add(TNamed('Description', "(Number of matching-hit layers) x 2 for true pions"))
451  ndof_pi.GetListOfFunctions().Add(TNamed('Check', "Even-only entries; sharp peak at 0 with tail above"))
452  ndof_pi.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
453  ndof_pi.GetListOfFunctions().Add(TNamed('MetaOptions', 'expert,pvalue-warn=0.50,pvalue-error=0.10'))
454  ndof_pi.GetXaxis().SetTitle('#dof')
455 
456  rchisq_mu.GetListOfFunctions().Add(TNamed('Description', "Transverse-deviation chi-squared/#dof for true muons"))
457  rchisq_mu.GetListOfFunctions().Add(TNamed('Check', "Narrow peak at 1"))
458  rchisq_mu.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
459  rchisq_mu.GetListOfFunctions().Add(TNamed('MetaOptions', 'expert,pvalue-warn=0.50,pvalue-error=0.10'))
460  rchisq_mu.GetXaxis().SetTitle('#chi^{2} / #ndof')
461  rchisq_pi.GetListOfFunctions().Add(TNamed('Description', "Transverse-deviation chi-squared/#dof for true pions"))
462  rchisq_pi.GetListOfFunctions().Add(TNamed('Check', "Broad peak at ~1"))
463  rchisq_pi.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
464  rchisq_pi.GetListOfFunctions().Add(TNamed('MetaOptions', 'expert,pvalue-warn=0.50,pvalue-error=0.10'))
465  rchisq_pi.GetXaxis().SetTitle('#chi^{2} / #ndof')
466 
467  efficiency_momentum.GetListOfFunctions().Add(TNamed('Description', "Muon ID efficiency vs p for logL(mu) > logL(pi)"))
468  efficiency_momentum.GetListOfFunctions().Add(TNamed('Check', "90-95% for p > 0.6 GeV/c"))
469  efficiency_momentum.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
470  efficiency_momentum.GetListOfFunctions().Add(TNamed('MetaOptions', 'pvalue-warn=0.50,pvalue-error=0.10'))
471  efficiency_momentum.GetXaxis().SetTitle('p (GeV/#it{c})')
472  efficiency_momentum.SetMinimum(0.0)
473  fakerate_momentum.GetListOfFunctions().Add(TNamed('Description', "Pion fake rate vs p for logL(mu) > logL(pi)"))
474  fakerate_momentum.GetListOfFunctions().Add(TNamed('Check', "Peak of 11% below 1 GeV/c, 2-6% above 1 GeV/c"))
475  fakerate_momentum.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
476  fakerate_momentum.GetListOfFunctions().Add(TNamed('MetaOptions', 'pvalue-warn=0.50,pvalue-error=0.10'))
477  fakerate_momentum.GetXaxis().SetTitle('p (GeV/#it{c})')
478  fakerate_momentum.SetMinimum(0.0)
479 
480  efficiency_theta.GetListOfFunctions().Add(TNamed('Description', "Muon ID efficiency vs theta for logL(mu) > logL(pi)"))
481  efficiency_theta.GetListOfFunctions().Add(TNamed('Check', "Somewhat flat (80-100%) in KLM range (17 < theta < 145 degrees)"))
482  efficiency_theta.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
483  efficiency_theta.GetListOfFunctions().Add(TNamed('MetaOptions', 'pvalue-warn=0.50,pvalue-error=0.10'))
484  efficiency_theta.GetXaxis().SetTitle('#theta (deg)')
485  efficiency_theta.SetMinimum(0.0)
486  fakerate_theta.GetListOfFunctions().Add(TNamed('Description', "Pion fake rate vs theta for logL(mu) > logL(pi)"))
487  fakerate_theta.GetListOfFunctions().Add(TNamed('Check', "Roughly flat (3-7%) in KLM range (17 < theta < 145 degrees)"))
488  fakerate_theta.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
489  fakerate_theta.GetListOfFunctions().Add(TNamed('MetaOptions', 'pvalue-warn=0.50,pvalue-error=0.10'))
490  fakerate_theta.GetXaxis().SetTitle('#theta (deg)')
491  fakerate_theta.SetMinimum(0.0)
492 
493  efficiency_phi.GetListOfFunctions().Add(TNamed('Description', "Muon identification efficiency vs phi for logL(mu) > logL(pi)"))
494  efficiency_phi.GetListOfFunctions().Add(TNamed('Check', "Somewhat flat distribution (85-100%)"))
495  efficiency_phi.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
496  efficiency_phi.GetListOfFunctions().Add(TNamed('MetaOptions', 'pvalue-warn=0.50,pvalue-error=0.10'))
497  efficiency_phi.GetXaxis().SetTitle('#phi (deg)')
498  efficiency_phi.SetMinimum(0.0)
499  fakerate_phi.GetListOfFunctions().Add(TNamed('Description', "Pion fake rate vs phi for logL(mu) > logL(pi)"))
500  fakerate_phi.GetListOfFunctions().Add(TNamed('Check', "Roughly flat distribution (2-7%)"))
501  fakerate_phi.GetListOfFunctions().Add(TNamed('Contact', contact_mail))
502  fakerate_phi.GetListOfFunctions().Add(TNamed('MetaOptions', 'pvalue-warn=0.50,pvalue-error=0.10'))
503  fakerate_phi.GetXaxis().SetTitle('#phi (deg)')
504  fakerate_phi.SetMinimum(0.0)
505 
506  outcome.Write()
507  llMu_mu.Write()
508  llMu_pi.Write()
509  llPi_mu.Write()
510  llPi_pi.Write()
511  llDiff_mu.Write()
512  llDiff_pi.Write()
513  barrelLayer_mu.Write()
514  barrelLayer_pi.Write()
515  endcapLayer_mu.Write()
516  endcapLayer_pi.Write()
517  diffLayer_mu.Write()
518  diffLayer_pi.Write()
519  ndof_mu.Write()
520  ndof_pi.Write()
521  rchisq_mu.Write()
522  rchisq_pi.Write()
523  efficiency_momentum.Write()
524  efficiency_theta.Write()
525  efficiency_phi.Write()
526  fakerate_momentum.Write()
527  fakerate_theta.Write()
528  fakerate_phi.Write()
529 
530 
531 
534 if __name__ == '__main__':
535  if ACTIVE:
536  main()
537  else:
538  print("This validation deactivated and thus basf2 is not executed.\n"
539  "If you want to run this validation, please set the 'ACTIVE' flag above to 'True'.\n"
540  "Exiting.")
Definition: main.py:1
int main(int argc, char **argv)
Run all tests.
Definition: test_main.cc:91