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