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>
35 ROOT.PyConfig.IgnoreCommandLineOptions =
True
36 from ROOT
import TFile, TChain, TH1F, \
37 gStyle, TNamed, TProfile
41 from optparse
import OptionParser
46 PT_VALUES = get_generated_pt_values()
51 CONTACT_PERSON = {
'Name':
'tracking software mailing list',
52 'Email':
'software-tracking@belle2.org'}
55 def make_expert_plot(hist):
56 hist.GetListOfFunctions().Add(TNamed(
'MetaOptions',
'expert'))
60 """Function which is executed"""
62 print(
'Tracking validation plots.')
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.'
69 option_parser.add_option(
'-o',
'--output-file', dest=
'output_file',
70 default=
'TrackingValidation.root',
71 help=
'Root file with all histograms.')
73 options = option_parser.parse_args()[0]
79 data_tree = TChain(tree_name)
80 data_tree.Add(options.input_file)
84 number_entries = data_tree.GetEntries()
85 except AttributeError:
86 print(f
'Could not load tree with name {tree_name}.')
88 if number_entries == 0:
89 print(f
'Data tree \'{tree_name}\'is empty or does not exist. Exit.')
93 output_file_name = options.output_file
94 output_root_file = TFile(output_file_name,
'recreate')
97 calculate_efficiency_in_pt(data_tree)
99 pt_of_interest = [0.05, 0.25, 1.]
102 for pt_value
in pt_of_interest:
104 generate_cos_theta_plot(data_tree, pt_value)
107 create_momentum_resolution_plot(data_tree)
109 draw_impact_parameter(data_tree)
111 draw_pvalue(data_tree)
113 draw_hit_counts(data_tree, pt_of_interest)
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)
127 pt_of_interest=pt_of_interest,
131 sub_dir =
'additionalPlots'
132 output_root_file.mkdir(sub_dir)
133 output_root_file.cd(sub_dir)
135 for pt_value
in PT_VALUES:
137 generate_cos_theta_plot(data_tree, pt_value)
139 draw_residua(data_tree,
'x',
'x_gen')
140 draw_residua(data_tree,
'y',
'y_gen')
141 draw_residua(data_tree,
'z',
'z_gen')
154 output_root_file.Write()
155 output_root_file.Close()
158 def draw_hit_counts(data_tree, pt_values):
160 Draw the hit count distribution.
166 number_bins = int((pt_upper - pt_lower) / 0.05 + 0.5)
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)
174 description = f
'Distribution of Hit Counts in {det} (Contact: {CONTACT_PERSON["Email"]}).'
177 hists[det].GetListOfFunctions().Add(TNamed(
'Description', description))
178 hists[det].GetListOfFunctions().Add(TNamed(
'Check', check))
179 make_expert_plot(hists[det])
181 hNweights = TH1F(
'hNweights',
'number of weights stored', 201, 0, 201)
182 hWeightsProfile = TProfile(
'hWeightsProfile',
'profile of weights', 201,
184 hWeightsProfile.SetMinimum(0)
185 hWeightsProfile.SetMaximum(1.1)
186 hWeights = TH1F(
'hWeights',
'DAF weights', 50, 0, 1)
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])
202 make_expert_plot(hWeightsProfile)
203 hWeightsProfile.Write()
206 def draw_pvalue(data_tree):
208 Create a histogram of the pValue of the track fits
211 number_entries = data_tree.GetEntries()
213 numBins = math.trunc(number_entries / 200)
214 hist_pvalue = TH1F(
'hpValue',
'p-Value Distribution', numBins, 0., 1.)
216 hist_pvalue.SetMinimum(0)
219 for entry
in range(number_entries):
220 data_tree.GetEntry(entry)
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.')
232 hist_pvalue.Fill(pvalue)
234 hist_pvalue.SetTitle(
'Distribution of pValue')
235 hist_pvalue.SetXTitle(
'pValue')
236 hist_pvalue.SetYTitle(
'number of entries')
238 description = f
'Distribution of pValues of the tracks (Contact: {CONTACT_PERSON["Email"]}).'
239 check =
'Should be a flat distribution.'
241 hist_pvalue.GetListOfFunctions().Add(TNamed(
'Description', description))
242 hist_pvalue.GetListOfFunctions().Add(TNamed(
'Check', check))
243 make_expert_plot(hist_pvalue)
248 def calculate_efficiency_in_pt(data_tree):
249 """Calculate single track reconstruction efficiency in bins of pt"""
254 number_bins = int((pt_upper - pt_lower) / 0.05 + 0.5)
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)
260 data_tree.Draw(
'pt_gen>>hPtGen',
'',
'goff')
261 data_tree.Draw(
'pt_gen>>hPtRec',
'pt != -999',
'goff')
264 efficiency_hist = TH1F(
'hEfficiency',
'hEfficiency', number_bins,
266 efficiency_hist.SetMinimum(0)
267 efficiency_hist.SetMaximum(1.05)
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.')
276 check =
'The efficiency should be stable for higher pt values.'
277 efficiency_hist.GetListOfFunctions().Add(TNamed(
'Description',
279 efficiency_hist.GetListOfFunctions().Add(TNamed(
'Check', check))
280 efficiency_hist.GetListOfFunctions().Add(TNamed(
'Contact',
281 CONTACT_PERSON[
'Email']))
284 for ibin
in range(1, number_bins + 1):
285 number_generated = hist_pt_gen.GetBinContent(ibin)
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))
294 efficiency_hist.SetBinContent(ibin, efficiency)
295 efficiency_hist.SetBinError(ibin, error)
297 efficiency_hist.SetTitle(
'Tracking efficiency in bins of transverse momentum pt.'
299 efficiency_hist.SetXTitle(
'pt in GeV')
300 efficiency_hist.SetYTitle(
'efficiency')
303 efficiency_hist.Write()
306 def generate_cos_theta_plot(data_tree, pt_value):
307 """Creates a efficiency histo in bins of cos theta"""
313 hist_cos_gen = TH1F(
'hCosGen',
'hCosGen', number_bins, cos_lower,
315 hist_cos_rec = TH1F(
'hCosRec',
'hCosRec', number_bins, cos_lower,
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')
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.'
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',
337 efficiency_hist.GetListOfFunctions().Add(TNamed(
'Check', check))
338 efficiency_hist.GetListOfFunctions().Add(TNamed(
'Contact',
339 CONTACT_PERSON[
'Email']))
341 for ibin
in range(1, number_bins + 1):
344 number_generated = hist_cos_gen.GetBinContent(ibin)
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))
351 efficiency_hist.SetBinContent(ibin, efficiency)
352 efficiency_hist.SetBinError(ibin, error)
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()
361 def create_momentum_resolution_plot(data_tree):
362 """Create momentum resolution plot"""
367 number_bins = int((pt_upper - pt_lower) / 0.05 + 0.5)
369 sigma_pt_values = calculate_momentum_resolution2(data_tree)
371 hist_resolution = TH1F(
'hPtResolution',
'hPtResolution', number_bins,
375 fit_pt_res_values = []
376 fit_pt_res_errors = []
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)
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
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)
391 hist_resolution.SetBinContent(ibin, sigma_pt_over_pt)
392 hist_resolution.SetBinError(ibin, sigma_pt_error_over_pt)
394 hist_resolution.SetTitle(
'Transverse Momentum resolution')
395 hist_resolution.SetXTitle(
'pt in GeV/c')
396 hist_resolution.SetYTitle(
'#sigma_{pt}/pt')
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.')
406 hist_resolution.GetListOfFunctions().Add(TNamed(
'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()
417 def get_scaled_rms_90(values, scale_factor=0.79):
421 param values: list of numbers
422 return scaled RMS90 of the distribution in values
425 (vMin, vMax) = (np.nanpercentile(values, 5), np.nanpercentile(values, 95))
426 final_list = np.array(values)[np.logical_and(values > vMin, values < vMax)]
428 rms_90 = np.sqrt(np.mean(np.square(final_list))) / scale_factor
430 rms_error = rms_90 / (scale_factor * np.sqrt(2 * len(final_list)))
432 return (rms_90, rms_error)
435 def value_in_list(test_value, value_list):
437 Function checks, if a given value is in the list values
440 for value
in value_list:
441 if value - DELTA_PT < test_value < value + DELTA_PT:
447 def calculate_momentum_resolution2(data_tree):
449 Calculate momentum resolution
454 for key
in PT_VALUES:
455 delta_pt_dict[key] = []
457 for entry_index
in range(data_tree.GetEntries()):
459 data_tree.GetEntry(entry_index)
461 pt_gen = data_tree.GetLeaf(
'pt_gen').GetValue()
463 test_object = value_in_list(pt_gen, PT_VALUES)
466 pt = data_tree.GetLeaf(
'pt').GetValue()
469 delta_pt_dict[test_object[1]].append(pt - pt_gen)
471 print(f
'pt = {test_object[1]:0.2f} is not generated and not used as key. Abort!')
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}')
482 return pt_resolutions
492 Fit a polynomial to data
494 @param x_value_list list of x values
496 @param y_value_list list of y values
498 @param degree degree of the fitting polynomial
500 @param y_value_errors list of corresponding errors on y x_values
504 import scipy.optimize
506 print(
'Module scipy.optimize is not installed. Fit cannot be done.')
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)
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)
525 print(par_opt, par_cov)
530 def draw_impact_parameter(data_tree):
532 Create a histogram with the impact parameter resolution
534 @param data_tree ROOT tree with the input data
537 print(
'********************************')
538 print(
'Impact parameter:')
542 for key
in PT_VALUES:
543 impact_param_d[key] = []
544 impact_param_z[key] = []
547 number_of_entries = data_tree.GetEntries()
548 except AttributeError:
549 print(
'ERROR: data_tree in function "draw_impact_parameter" is no TTree.')
552 for index_entry
in range(number_of_entries):
553 data_tree.GetEntry(index_entry)
555 pt_gen = data_tree.GetLeaf(
'pt_gen').GetValue()
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()
564 poca_y = data_tree.GetLeaf(
'y').GetValue()
565 poca_z = data_tree.GetLeaf(
'z').GetValue()
567 z_gen = data_tree.GetLeaf(
'z_gen').GetValue()
569 px = data_tree.GetLeaf(
'px').GetValue()
570 py = data_tree.GetLeaf(
'py').GetValue()
572 d0 = d0_sign([poca_x, poca_y], [px, py]) * np.hypot(poca_x, poca_y)
574 impact_param_d[test_object[1]].append(d0)
575 impact_param_z[test_object[1]].append(poca_z - z_gen)
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}')
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']))
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']))
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])
618 hist_impact_parameter_d0.Write()
619 hist_impact_parameter_z.Write()
622 def scale_histogram(hist):
624 Scales a TH1F using the TH1::Scale function
627 if isinstance(hist, TH1F):
628 if hist.GetEntries() > 0:
629 hist.Scale(1 / hist.GetEntries())
631 print(
'Cannot scale. No entries in histogram ' + hist.GetName())
633 print(
'Not a TH1F. Not able to scale.')
636 def d0_sign(vector_d, momentum):
638 Estimate the sign of the impact parameter d0_sign
643 if len(vector_d) != 2
or len(momentum) != 2:
644 print(
'Need x and y values of d and p_d.')
647 phi_momentum = np.arctan2(momentum[1], momentum[0])
648 phi_d = np.arctan2(vector_d[1], vector_d[0])
650 diff = phi_momentum - phi_d
653 diff = diff - 2 * np.pi
655 diff = 2 * np.pi + diff
671 Write histogram of difference of variable_name (reconstructed) and
672 gen_variable_name (generated) in gDirectory
675 if pt_of_interest
is None:
678 used_pts = [pt
for pt
in pt_of_interest
if pt
in PT_VALUES]
680 number_entries = data_tree.GetEntries()
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',
688 for index
in range(number_entries):
689 data_tree.GetEntry(index)
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:
698 histograms[pt_gen].Fill((variable_reconstructed - variable_generated) / variable_generated)
700 histograms[pt_gen].Fill(variable_reconstructed - variable_generated)
702 for (pt, hist)
in histograms.items():
703 scale_histogram(hist)
705 hist.SetXTitle(f
'({variable_name} - {gen_variable_name}) / {gen_variable_name}')
707 hist.SetXTitle(f
'({variable_name} - {gen_variable_name})')
708 hist.GetListOfFunctions().Add(TNamed(
'Contact', CONTACT_PERSON[
'Email'
710 make_expert_plot(hist)
714 def tf1_fit_function(x, par):
715 return fit_function_sigma_pt_over_pt(x, par[0], par[1])
718 def fit_function_sigma_pt_over_pt(x, a, b):
719 return np.sqrt((a * x) ** 2 + b * b)
722 if __name__ ==
'__main__':
int main(int argc, char **argv)
Run all tests.