Belle II Software  release-08-00-04
13_trackingEfficiency_createPlots.py
1 #!/usr/bin/env python3
2 
3 
10 
11 
17 
18 """
19 <header>
20  <input>trackingEfficiency_pt_0.05GeV.root,
21  trackingEfficiency_pt_0.10GeV.root,trackingEfficiency_pt_0.25GeV.root,
22  trackingEfficiency_pt_0.40GeV.root,trackingEfficiency_pt_0.60GeV.root,
23  trackingEfficiency_pt_1.00GeV.root,trackingEfficiency_pt_1.50GeV.root,
24  trackingEfficiency_pt_2.00GeV.root,trackingEfficiency_pt_3.00GeV.root,
25  trackingEfficiency_pt_4.00GeV.root</input>
26  <output>TrackingValidation.root</output>
27  <contact>software-tracking@belle2.org</contact>
28  <description>Create momentum resolution, impact parameter resolution and efficiency plots.</description>
29 </header>
30 """
31 
32 
33 import ROOT
34 
35 ROOT.PyConfig.IgnoreCommandLineOptions = True # noqa
36 from ROOT import TFile, TChain, TH1F, \
37  gStyle, TNamed, TProfile
38 import sys
39 import math
40 import numpy as np
41 from optparse import OptionParser
42 from tracking.validation.tracking_efficiency_helpers import get_generated_pt_values
43 
44 DELTA_PT = 0.0001
45 
46 PT_VALUES = get_generated_pt_values()
47 # PT_VALUES = [0.05, 0.1, 0.15, 0.25, 0.5, 1.0, 2.0, 3.0]
48 
49 # contact person information
50 # is added to the plot descriptions
51 CONTACT_PERSON = {'Name': 'tracking software mailing list',
52  'Email': 'software-tracking@belle2.org'}
53 
54 
55 def make_expert_plot(hist):
56  hist.GetListOfFunctions().Add(TNamed('MetaOptions', 'expert'))
57 
58 
59 def main():
60  """Function which is executed"""
61 
62  print('Tracking validation plots.')
63 
64  option_parser = OptionParser()
65  option_parser.add_option('-i', '--input-file', dest='input_file',
66  default='../trackingEfficiency_pt_*.root',
67  help='Root file with StandardTrackingPerformance output.'
68  )
69  option_parser.add_option('-o', '--output-file', dest='output_file',
70  default='TrackingValidation.root',
71  help='Root file with all histograms.')
72 
73  options = option_parser.parse_args()[0]
74 
75  gStyle.SetOptStat(0)
76 
77  # load data tree
78  tree_name = 'data'
79  data_tree = TChain(tree_name) # input_root_file.Get(tree_name)
80  data_tree.Add(options.input_file)
81 
82  number_entries = 0
83  try:
84  number_entries = data_tree.GetEntries()
85  except AttributeError:
86  print(f'Could not load tree with name {tree_name}.')
87 
88  if number_entries == 0:
89  print(f'Data tree \'{tree_name}\'is empty or does not exist. Exit.')
90  sys.exit(0)
91 
92  # output root file
93  output_file_name = options.output_file
94  output_root_file = TFile(output_file_name, 'recreate')
95 
96  # create efficiency in bins of pt plot
97  calculate_efficiency_in_pt(data_tree)
98 
99  pt_of_interest = [0.05, 0.25, 1.]
100  # pt_of_interest = PT_VALUES
101 
102  for pt_value in pt_of_interest:
103  # create plots of efficiency in bins of cos Theta for different pt
104  generate_cos_theta_plot(data_tree, pt_value)
105 
106  # create momentum resolution plot
107  create_momentum_resolution_plot(data_tree)
108 
109  draw_impact_parameter(data_tree)
110 
111  draw_pvalue(data_tree)
112 
113  draw_hit_counts(data_tree, pt_of_interest)
114 
115  draw_residua(data_tree, 'x', 'x_gen', pt_of_interest=pt_of_interest)
116  draw_residua(data_tree, 'y', 'y_gen', pt_of_interest=pt_of_interest)
117  draw_residua(data_tree, 'z', 'z_gen', pt_of_interest=pt_of_interest)
118 
119  draw_residua(
120  data_tree,
121  'pt',
122  'pt_gen',
123  bins=200,
124  ledge=-0.1,
125  uedge=0.1,
126  normalize=1,
127  pt_of_interest=pt_of_interest,
128  )
129 
130  # write additional plots in sub dir
131  sub_dir = 'additionalPlots'
132  output_root_file.mkdir(sub_dir)
133  output_root_file.cd(sub_dir)
134 
135  for pt_value in PT_VALUES:
136  # create plots of efficiency in bins of cos Theta for different pt
137  generate_cos_theta_plot(data_tree, pt_value)
138 
139  draw_residua(data_tree, 'x', 'x_gen')
140  draw_residua(data_tree, 'y', 'y_gen')
141  draw_residua(data_tree, 'z', 'z_gen')
142 
143  draw_residua(
144  data_tree,
145  'pt',
146  'pt_gen',
147  bins=200,
148  ledge=-0.1,
149  uedge=0.1,
150  normalize=1,
151  )
152 
153  # close output file
154  output_root_file.Write()
155  output_root_file.Close()
156 
157 
158 def draw_hit_counts(data_tree, pt_values):
159  """
160  Draw the hit count distribution.
161  """
162 
163  pt_lower = -0.025
164  pt_upper = 4.025
165  # 50 MeV wide bins. + 0.5 serves for correct rounding.
166  number_bins = int((pt_upper - pt_lower) / 0.05 + 0.5)
167 
168  hists = {}
169  for det in ['PXD', 'SVD', 'CDC']:
170  hists[det] = TProfile(f'hHitCounts{det}',
171  f'Hit count profile for the {det};pT;nHits',
172  number_bins, pt_lower, pt_upper)
173 
174  description = f'Distribution of Hit Counts in {det} (Contact: {CONTACT_PERSON["Email"]}).'
175  check = ''
176 
177  hists[det].GetListOfFunctions().Add(TNamed('Description', description))
178  hists[det].GetListOfFunctions().Add(TNamed('Check', check))
179  make_expert_plot(hists[det])
180 
181  hNweights = TH1F('hNweights', 'number of weights stored', 201, 0, 201)
182  hWeightsProfile = TProfile('hWeightsProfile', 'profile of weights', 201,
183  0, 201)
184  hWeightsProfile.SetMinimum(0)
185  hWeightsProfile.SetMaximum(1.1)
186  hWeights = TH1F('hWeights', 'DAF weights', 50, 0, 1)
187 
188  for track in data_tree:
189  if not track.pValue == -999:
190  hists['PXD'].Fill(track.pt_gen, track.nPXDhits)
191  hists['SVD'].Fill(track.pt_gen, track.nSVDhits)
192  hists['CDC'].Fill(track.pt_gen, track.nCDChits)
193  hNweights.Fill(track.nWeights)
194  for i in range(track.nWeights):
195  hWeights.Fill(track.weights[i])
196  hWeightsProfile.Fill(i, track.weights[i])
197 
198  for key in hists:
199  hists[key].Write()
200  # hNweights.Write()
201  # hWeights.Write()
202  make_expert_plot(hWeightsProfile)
203  hWeightsProfile.Write()
204 
205 
206 def draw_pvalue(data_tree):
207  """
208  Create a histogram of the pValue of the track fits
209  """
210 
211  number_entries = data_tree.GetEntries()
212  # Normalize number of bins to number of events to keep plots comparable.
213  numBins = math.trunc(number_entries / 200)
214  hist_pvalue = TH1F('hpValue', 'p-Value Distribution', numBins, 0., 1.)
215  # Make axis range independent of first bin contents for more pleasant look
216  hist_pvalue.SetMinimum(0)
217  # hist_pvalue.SetMaximum(number_entries / 50.)
218 
219  for entry in range(number_entries):
220  data_tree.GetEntry(entry)
221 
222  try:
223  pvalue = data_tree.GetLeaf('pValue').GetValue()
224  except ReferenceError:
225  print(f'The variable "pValue" doesn\'t exist in the tree "{data_tree.GetName()}".\n'
226  'Leave this function without plotting the variable.')
227  return
228 
229  if pvalue is -999:
230  continue
231  else:
232  hist_pvalue.Fill(pvalue)
233 
234  hist_pvalue.SetTitle('Distribution of pValue')
235  hist_pvalue.SetXTitle('pValue')
236  hist_pvalue.SetYTitle('number of entries')
237 
238  description = f'Distribution of pValues of the tracks (Contact: {CONTACT_PERSON["Email"]}).'
239  check = 'Should be a flat distribution.'
240 
241  hist_pvalue.GetListOfFunctions().Add(TNamed('Description', description))
242  hist_pvalue.GetListOfFunctions().Add(TNamed('Check', check))
243  make_expert_plot(hist_pvalue)
244 
245  hist_pvalue.Write()
246 
247 
248 def calculate_efficiency_in_pt(data_tree):
249  """Calculate single track reconstruction efficiency in bins of pt"""
250 
251  pt_lower = -0.025
252  pt_upper = 4.025
253  # 50 MeV wide bins. + 0.5 serves for correct rounding.
254  number_bins = int((pt_upper - pt_lower) / 0.05 + 0.5)
255 
256  hist_pt_gen = TH1F('hPtGen', 'hPtGen', number_bins, pt_lower, pt_upper)
257  hist_pt_rec = TH1F('hPtRec', 'hPtRec', number_bins, pt_lower, pt_upper)
258 
259  # draw data in defined histos
260  data_tree.Draw('pt_gen>>hPtGen', '', 'goff')
261  data_tree.Draw('pt_gen>>hPtRec', 'pt != -999', 'goff')
262 
263  # final hist
264  efficiency_hist = TH1F('hEfficiency', 'hEfficiency', number_bins,
265  pt_lower, pt_upper)
266  efficiency_hist.SetMinimum(0)
267  efficiency_hist.SetMaximum(1.05)
268 
269  description = ('Events with 10 muon tracks with fixed transverse '
270  'momentum are generated using the ParticleGun(500 '
271  'events for each pt value). The events are reconstructed '
272  'with VXDTF + CDC Track Finder + VXDCDCMerger.'
273  'plot shows the single track reconstruction efficiency '
274  'over the transverse momentum.')
275 
276  check = 'The efficiency should be stable for higher pt values.'
277  efficiency_hist.GetListOfFunctions().Add(TNamed('Description',
278  description))
279  efficiency_hist.GetListOfFunctions().Add(TNamed('Check', check))
280  efficiency_hist.GetListOfFunctions().Add(TNamed('Contact',
281  CONTACT_PERSON['Email']))
282 
283  # loop over bins and calculate efficiency and error of efficiency
284  for ibin in range(1, number_bins + 1):
285  number_generated = hist_pt_gen.GetBinContent(ibin)
286  efficiency = 0
287  error = 0
288 
289  if number_generated > 0:
290  number_reconstructed = hist_pt_rec.GetBinContent(ibin)
291  efficiency = number_reconstructed / number_generated
292  error = math.sqrt(number_reconstructed * (number_generated - number_reconstructed) / pow(number_generated, 3))
293  # set according bin to the efficiency value and add errorbar
294  efficiency_hist.SetBinContent(ibin, efficiency)
295  efficiency_hist.SetBinError(ibin, error)
296 
297  efficiency_hist.SetTitle('Tracking efficiency in bins of transverse momentum pt.'
298  )
299  efficiency_hist.SetXTitle('pt in GeV')
300  efficiency_hist.SetYTitle('efficiency')
301 
302  # write hist to the output file
303  efficiency_hist.Write()
304 
305 
306 def generate_cos_theta_plot(data_tree, pt_value):
307  """Creates a efficiency histo in bins of cos theta"""
308 
309  number_bins = 100
310  cos_lower = -1
311  cos_upper = 1
312 
313  hist_cos_gen = TH1F('hCosGen', 'hCosGen', number_bins, cos_lower,
314  cos_upper)
315  hist_cos_rec = TH1F('hCosRec', 'hCosRec', number_bins, cos_lower,
316  cos_upper)
317 
318  data_tree.Draw('cosTheta_gen>>hCosGen',
319  f'pt_gen>({pt_value:.2f} - {DELTA_PT:f}) && pt_gen<({pt_value:.2f} + {DELTA_PT:f})', 'goff')
320  data_tree.Draw('cosTheta_gen>>hCosRec',
321  f'pt_gen>({pt_value:.2f} - {DELTA_PT:f}) && pt_gen<({pt_value:.2f} + {DELTA_PT:f}) && pt != -999', 'goff')
322 
323  description = ('Events with 10 muon tracks with fixed transverse '
324  'momentum are generated using the ParticleGun(500 '
325  'events for each pt value). The events are reconstructed '
326  'with VXDTF + CDCTF + MCTrackCandCombiner. The plot '
327  'shows the single track reconstruction efficiency in '
328  'bins of the polar angle for the fixed transverse '
329  'momentum pt = %.2f GeV.')
330  check = 'Stable efficiency over the hole range of the polar angle.'
331 
332  efficiency_hist = TH1F(f'hEfficiencyPt{pt_value:.2f}GeV',
333  f'hEfficiencyPt{pt_value:.2f}GeV',
334  number_bins, cos_lower, cos_upper)
335  efficiency_hist.GetListOfFunctions().Add(TNamed('Description',
336  description))
337  efficiency_hist.GetListOfFunctions().Add(TNamed('Check', check))
338  efficiency_hist.GetListOfFunctions().Add(TNamed('Contact',
339  CONTACT_PERSON['Email']))
340 
341  for ibin in range(1, number_bins + 1):
342  efficiency = 0
343  error = 0
344  number_generated = hist_cos_gen.GetBinContent(ibin)
345 
346  if number_generated > 0:
347  number_reconstructed = hist_cos_rec.GetBinContent(ibin)
348  efficiency = number_reconstructed / number_generated
349  error = math.sqrt(number_reconstructed * (number_generated - number_reconstructed) / pow(number_generated, 3))
350 
351  efficiency_hist.SetBinContent(ibin, efficiency)
352  efficiency_hist.SetBinError(ibin, error)
353 
354  efficiency_hist.SetTitle(f'Tracks with pt = {pt_value:.2f} GeV')
355  efficiency_hist.SetXTitle('cos #Theta')
356  efficiency_hist.SetYTitle('efficiency')
357  make_expert_plot(efficiency_hist)
358  efficiency_hist.Write()
359 
360 
361 def create_momentum_resolution_plot(data_tree):
362  """Create momentum resolution plot"""
363 
364  pt_lower = -0.025
365  pt_upper = 4.025
366  # 50 MeV wide bins. + 0.5 serves for correct rounding.
367  number_bins = int((pt_upper - pt_lower) / 0.05 + 0.5)
368 
369  sigma_pt_values = calculate_momentum_resolution2(data_tree)
370 
371  hist_resolution = TH1F('hPtResolution', 'hPtResolution', number_bins,
372  pt_lower, pt_upper)
373 
374  fit_pt_values = []
375  fit_pt_res_values = []
376  fit_pt_res_errors = []
377 
378  for (pt_value, sigma_pt_value) in sigma_pt_values.items():
379  ibin = hist_resolution.FindBin(pt_value)
380  bin_center = hist_resolution.GetBinCenter(ibin)
381 
382  sigma_pt_over_pt = 0
383  if sigma_pt_value[0] != -999:
384  sigma_pt_over_pt = sigma_pt_value[0] / pt_value
385  sigma_pt_error_over_pt = sigma_pt_value[1] / pt_value
386 
387  fit_pt_res_values.append(sigma_pt_over_pt)
388  fit_pt_res_errors.append(sigma_pt_error_over_pt)
389  fit_pt_values.append(bin_center)
390 
391  hist_resolution.SetBinContent(ibin, sigma_pt_over_pt)
392  hist_resolution.SetBinError(ibin, sigma_pt_error_over_pt)
393 
394  hist_resolution.SetTitle('Transverse Momentum resolution')
395  hist_resolution.SetXTitle('pt in GeV/c')
396  hist_resolution.SetYTitle('#sigma_{pt}/pt')
397 
398  description = ('Events with 10 muon tracks with fixed transverse '
399  'momentum are generated using the ParticleGun(200 '
400  'events for each pt value). The events are reconstructed '
401  'with VXDTF + CDCTF + MCTrackCandCombiner. The plot '
402  'shows the relative momentum resolution of the '
403  'transverse momentum over transverse momentum.')
404 
405  check = ''
406  hist_resolution.GetListOfFunctions().Add(TNamed('Description',
407  description))
408  hist_resolution.GetListOfFunctions().Add(TNamed('Check', check))
409  hist_resolution.GetListOfFunctions().Add(TNamed('Contact',
410  CONTACT_PERSON['Email']))
411  hist_resolution.GetListOfFunctions().Add(TNamed('MetaOptions', 'logy'))
412  hist_resolution.Write()
413 
414  # fit_parameters = fit_resolution(fit_pt_values, fit_pt_res_values, 1, fit_pt_res_errors)
415 
416 
417 def get_scaled_rms_90(values, scale_factor=0.79):
418  """
419  Calculate the RMS90
420 
421  param values: list of numbers
422  return scaled RMS90 of the distribution in values
423  """
424 
425  (vMin, vMax) = (np.nanpercentile(values, 5), np.nanpercentile(values, 95))
426  final_list = np.array(values)[np.logical_and(values > vMin, values < vMax)]
427 
428  rms_90 = np.sqrt(np.mean(np.square(final_list))) / scale_factor
429 
430  rms_error = rms_90 / (scale_factor * np.sqrt(2 * len(final_list)))
431 
432  return (rms_90, rms_error)
433 
434 
435 def value_in_list(test_value, value_list):
436  """
437  Function checks, if a given value is in the list values
438  """
439 
440  for value in value_list:
441  if value - DELTA_PT < test_value < value + DELTA_PT:
442  return (True, value)
443 
444  return (False, -999)
445 
446 
447 def calculate_momentum_resolution2(data_tree):
448  """
449  Calculate momentum resolution
450  """
451 
452  delta_pt_dict = {}
453 
454  for key in PT_VALUES:
455  delta_pt_dict[key] = []
456 
457  for entry_index in range(data_tree.GetEntries()):
458  # load entry
459  data_tree.GetEntry(entry_index)
460 
461  pt_gen = data_tree.GetLeaf('pt_gen').GetValue()
462 
463  test_object = value_in_list(pt_gen, PT_VALUES)
464 
465  if test_object[0]:
466  pt = data_tree.GetLeaf('pt').GetValue()
467  if pt != -999:
468  try:
469  delta_pt_dict[test_object[1]].append(pt - pt_gen)
470  except KeyError:
471  print(f'pt = {test_object[1]:0.2f} is not generated and not used as key. Abort!')
472  sys.exit(1)
473  pt_resolutions = {}
474 
475  print('********************************')
476  print('Momentum resolution:')
477  for key in sorted(delta_pt_dict):
478  (rms90, rms_error) = get_scaled_rms_90(delta_pt_dict[key])
479  pt_resolutions[key] = [rms90, rms_error]
480  print(f'pt = {key:0.2f}: sigma/pt = {rms90/key:0.4f}')
481 
482  return pt_resolutions
483 
484 
485 def fit_resolution(
486  x_value_list,
487  y_value_list,
488  degree,
489  y_value_errors=None,
490 ):
491  """
492  Fit a polynomial to data
493 
494  @param x_value_list list of x values
495 
496  @param y_value_list list of y values
497 
498  @param degree degree of the fitting polynomial
499 
500  @param y_value_errors list of corresponding errors on y x_values
501  """
502 
503  try:
504  import scipy.optimize
505  except ImportError:
506  print('Module scipy.optimize is not installed. Fit cannot be done.')
507  return []
508 
509  # weights = None
510  # if y_value_errors is not None:
511  # weights = [1 / x ** 2 for x in y_value_errors]
512  #
513  # parameters = np.polyfit(x_value_list, y_value_list, degree, w=weights,
514  # cov=True)
515 
516  x_value_list = np.array(x_value_list)
517  y_value_list = np.array(y_value_list)
518  y_value_errors = np.array(y_value_errors)
519 
520  (par_opt, par_cov) = \
521  scipy.optimize.curve_fit(fit_function_sigma_pt_over_pt, x_value_list,
522  y_value_list, [0.002, 0.003],
523  sigma=y_value_errors)
524 
525  print(par_opt, par_cov)
526 
527  return par_opt
528 
529 
530 def draw_impact_parameter(data_tree):
531  """
532  Create a histogram with the impact parameter resolution
533 
534  @param data_tree ROOT tree with the input data
535  """
536 
537  print('********************************')
538  print('Impact parameter:')
539 
540  impact_param_d = {}
541  impact_param_z = {}
542  for key in PT_VALUES:
543  impact_param_d[key] = []
544  impact_param_z[key] = []
545 
546  try:
547  number_of_entries = data_tree.GetEntries()
548  except AttributeError:
549  print('ERROR: data_tree in function "draw_impact_parameter" is no TTree.')
550  sys.exit(2)
551 
552  for index_entry in range(number_of_entries):
553  data_tree.GetEntry(index_entry)
554 
555  pt_gen = data_tree.GetLeaf('pt_gen').GetValue()
556 
557  test_object = value_in_list(pt_gen, PT_VALUES)
558  if test_object[0] is True:
559  poca_x = data_tree.GetLeaf('x').GetValue()
560 
561  if poca_x == -999:
562  continue
563 
564  poca_y = data_tree.GetLeaf('y').GetValue()
565  poca_z = data_tree.GetLeaf('z').GetValue()
566 
567  z_gen = data_tree.GetLeaf('z_gen').GetValue()
568 
569  px = data_tree.GetLeaf('px').GetValue()
570  py = data_tree.GetLeaf('py').GetValue()
571 
572  d0 = d0_sign([poca_x, poca_y], [px, py]) * np.hypot(poca_x, poca_y)
573 
574  impact_param_d[test_object[1]].append(d0)
575  impact_param_z[test_object[1]].append(poca_z - z_gen)
576  # estimate sign of d0
577 
578  d0_resolutions = {}
579  z_resolutions = {}
580  for key in sorted(impact_param_d):
581  (d0_resolution, d0_error) = get_scaled_rms_90(impact_param_d[key])
582  d0_resolutions[key] = [d0_resolution, d0_error]
583  (z_resolution, z_error) = get_scaled_rms_90(impact_param_z[key])
584  z_resolutions[key] = [z_resolution, z_error]
585  print(f'pt = {key:0.2f}: sigma_d0 = {d0_resolution:0.4f}, sigma_z = {z_resolution:0.4f}')
586 
587  number_bins = 62
588  lower_edge = -0.025
589  upper_edge = 4.025
590 
591  hist_impact_parameter_d0 = TH1F('hd0resolution', 'hd0resolution',
592  number_bins, lower_edge, upper_edge)
593  hist_impact_parameter_d0.SetTitle('d0 resolution')
594  hist_impact_parameter_d0.SetXTitle('pt in GeV/c')
595  hist_impact_parameter_d0.SetYTitle('#sigma_{d0} [cm]')
596  hist_impact_parameter_d0.GetListOfFunctions().Add(TNamed('MetaOptions', 'logy'))
597  hist_impact_parameter_d0.GetListOfFunctions().Add(TNamed('Contact',
598  CONTACT_PERSON['Email']))
599 
600  hist_impact_parameter_z = TH1F('hzresolution', 'hzresolution',
601  number_bins, lower_edge, upper_edge)
602  hist_impact_parameter_z.SetTitle('z resolution')
603  hist_impact_parameter_z.SetXTitle('pt in GeV/c')
604  hist_impact_parameter_z.SetYTitle('#sigma_{z} [cm]')
605  hist_impact_parameter_z.GetListOfFunctions().Add(TNamed('MetaOptions', 'logy'))
606  hist_impact_parameter_z.GetListOfFunctions().Add(TNamed('Contact',
607  CONTACT_PERSON['Email']))
608 
609  for pt in PT_VALUES:
610  ibin = hist_impact_parameter_d0.FindBin(pt)
611  if d0_resolutions[pt][0] > 0 and d0_resolutions[pt][1] > 0:
612  hist_impact_parameter_d0.SetBinContent(ibin, d0_resolutions[pt][0])
613  hist_impact_parameter_d0.SetBinError(ibin, d0_resolutions[pt][1])
614  if z_resolutions[pt][0] > 0 and z_resolutions[pt][1] > 0:
615  hist_impact_parameter_z.SetBinContent(ibin, z_resolutions[pt][0])
616  hist_impact_parameter_z.SetBinError(ibin, z_resolutions[pt][1])
617 
618  hist_impact_parameter_d0.Write()
619  hist_impact_parameter_z.Write()
620 
621 
622 def scale_histogram(hist):
623  """
624  Scales a TH1F using the TH1::Scale function
625  """
626 
627  if isinstance(hist, TH1F):
628  if hist.GetEntries() > 0:
629  hist.Scale(1 / hist.GetEntries())
630  else:
631  print('Cannot scale. No entries in histogram ' + hist.GetName())
632  else:
633  print('Not a TH1F. Not able to scale.')
634 
635 
636 def d0_sign(vector_d, momentum):
637  """
638  Estimate the sign of the impact parameter d0_sign
639 
640  @return -1 or +1
641  """
642 
643  if len(vector_d) != 2 or len(momentum) != 2:
644  print('Need x and y values of d and p_d.')
645  sys.exit(1)
646 
647  phi_momentum = np.arctan2(momentum[1], momentum[0])
648  phi_d = np.arctan2(vector_d[1], vector_d[0])
649 
650  diff = phi_momentum - phi_d
651  # Catch corner cases. For clarity, don't return early.
652  if diff > np.pi:
653  diff = diff - 2 * np.pi
654  elif diff < -np.pi:
655  diff = 2 * np.pi + diff
656 
657  return np.sign(diff)
658 
659 
660 def draw_residua(
661  data_tree,
662  variable_name,
663  gen_variable_name,
664  bins=100,
665  ledge=-0.01,
666  uedge=0.01,
667  normalize=0,
668  pt_of_interest=None,
669 ):
670  """
671  Write histogram of difference of variable_name (reconstructed) and
672  gen_variable_name (generated) in gDirectory
673  """
674 
675  if pt_of_interest is None:
676  used_pts = PT_VALUES
677  else:
678  used_pts = [pt for pt in pt_of_interest if pt in PT_VALUES]
679 
680  number_entries = data_tree.GetEntries()
681 
682  histograms = {}
683  for pt_value in used_pts:
684  histograms[pt_value] = TH1F(f'h{variable_name}Residuum_{pt_value:0.2f}GeV',
685  f'h{variable_name}Residuum_{pt_value:0.2f}GeV',
686  bins, ledge, uedge)
687 
688  for index in range(number_entries):
689  data_tree.GetEntry(index)
690 
691  variable_reconstructed = data_tree.GetLeaf(variable_name).GetValue()
692  if variable_reconstructed != -999:
693  variable_generated = \
694  data_tree.GetLeaf(gen_variable_name).GetValue()
695  pt_gen = round(data_tree.GetLeaf('pt_gen').GetValue(), 2)
696  if pt_gen in used_pts:
697  if normalize:
698  histograms[pt_gen].Fill((variable_reconstructed - variable_generated) / variable_generated)
699  else:
700  histograms[pt_gen].Fill(variable_reconstructed - variable_generated)
701 
702  for (pt, hist) in histograms.items():
703  scale_histogram(hist)
704  if normalize:
705  hist.SetXTitle(f'({variable_name} - {gen_variable_name}) / {gen_variable_name}')
706  else:
707  hist.SetXTitle(f'({variable_name} - {gen_variable_name})')
708  hist.GetListOfFunctions().Add(TNamed('Contact', CONTACT_PERSON['Email'
709  ]))
710  make_expert_plot(hist)
711  hist.Write()
712 
713 
714 def tf1_fit_function(x, par):
715  return fit_function_sigma_pt_over_pt(x, par[0], par[1])
716 
717 
718 def fit_function_sigma_pt_over_pt(x, a, b):
719  return np.sqrt((a * x) ** 2 + b * b)
720 
721 
722 if __name__ == '__main__':
723  main()
Definition: main.py:1
int main(int argc, char **argv)
Run all tests.
Definition: test_main.cc:91