Belle II Software  release-06-01-15
v0ValidationCreatePlots.py
1 #!/usr/bin/env python3
2 # -*- coding: utf-8 -*-
3 
4 
11 
12 """
13 <header>
14  <contact>software-tracking@belle2.org</contact>
15  <input>V0ValidationHarvested.root</input>
16  <description>This module creates efficiency plots for the V0 validation.</description>
17 </header>
18 """
19 
20 
21 import numpy
22 import ROOT
23 
24 
26 
27  """Reads the output created by the V0Harvester and creates plots from it."""
28 
29  def __init__(self, input_file='../V0ValidationHarvested.root', output_file='V0Validation.root'):
30  """Reads the output created by the V0Harvester defines histograms which will be filled later.
31 
32  :param input_file: Output of V0ValidationHarvester.
33  :param output_file: Plots displayed in the V0Validation.
34  """
35 
36  self.input_fileinput_file = input_file
37 
38  self.output_fileoutput_file = output_file
39 
40 
41  self.hist_rhist_r = ROOT.TH1F("", "True R", 20, 0, 20)
42 
43  self.hist_thetahist_theta = ROOT.TH1F("", "True Theta", 26, 20, 150)
44 
45  self.hist_phihist_phi = ROOT.TH1F("", "True Phi", 36, -180, 180)
46 
47  self.hist_phist_p = ROOT.TH1F("", "True P", 25, 0.0, 1.0)
48 
49 
50  self.hist_r_foundhist_r_found = ROOT.TH1F("", "Found R", 20, 0, 20)
51 
52  self.hist_theta_foundhist_theta_found = ROOT.TH1F("", "Found Theta", 26, 20, 150)
53 
54  self.hist_phi_foundhist_phi_found = ROOT.TH1F("", "Found Phi", 36, -180, 180)
55 
56  self.hist_p_foundhist_p_found = ROOT.TH1F("", "Found P", 25, 0.0, 1.0)
57 
58 
59  self.hist_invariant_masshist_invariant_mass = ROOT.TH1F("", "", 60, 0.470, 0.530)
60 
61  self.hist_invariant_mass_reshist_invariant_mass_res = ROOT.TH1F("", "", 40, -0.02, 0.02)
62 
63 
64  self.hist_chi2hist_chi2 = ROOT.TH1F("", "", 50, 0, 50)
65 
66  self.hist_chi2_insidehist_chi2_inside = ROOT.TH1F("", "", 50, 0, 50)
67 
68  self.hist_chi2_outsidehist_chi2_outside = ROOT.TH1F("", "", 50, 0, 50)
69 
70 
71  self.hist_mass_vs_mc_masshist_mass_vs_mc_mass = ROOT.TH2F("", "", 80, 0, 0.8, 80, 0, 0.8)
72 
73 
74  self.hist_invariant_mass_residuumhist_invariant_mass_residuum = ROOT.TH1F("", "", 60, -0.05, 0.05)
75 
76  self.hist_r_residuumhist_r_residuum = ROOT.TH1F("", "", 60, -0.1, 0.1)
77 
78  self.hist_theta_residuumhist_theta_residuum = ROOT.TH1F("", "", 60, -0.1, 0.1)
79 
80  self.hist_phi_residuumhist_phi_residuum = ROOT.TH1F("", "", 60, -0.1, 0.1)
81 
82  self.hist_p_residuumhist_p_residuum = ROOT.TH1F("", "", 60, -0.05, 0.05)
83 
84  def collect_histograms(self):
85  """Fills the defined histograms with the V0Harvester data.
86 
87  :return: self
88  """
89  input_root_file = ROOT.TFile.Open(self.input_fileinput_file, "READ")
90 
91  for event in input_root_file.V0Harvester_tree:
92  self.hist_rhist_r.Fill(event.R_MC)
93  self.hist_thetahist_theta.Fill(numpy.rad2deg(event.THETA_MC))
94  self.hist_phihist_phi.Fill(numpy.rad2deg(event.PHI_MC))
95  self.hist_phist_p.Fill(event.P_MC)
96 
97  if event.FOUND:
98  self.hist_r_foundhist_r_found.Fill(event.R_MC)
99  self.hist_theta_foundhist_theta_found.Fill(numpy.rad2deg(event.THETA_MC))
100  self.hist_phi_foundhist_phi_found.Fill(numpy.rad2deg(event.PHI_MC))
101  self.hist_p_foundhist_p_found.Fill(event.P_MC)
102 
103  self.hist_invariant_masshist_invariant_mass.Fill(event.M)
104 
105  # Residuum histograms
106  self.hist_invariant_mass_residuumhist_invariant_mass_residuum.Fill(event.M - event.M_MC)
107  self.hist_r_residuumhist_r_residuum.Fill(event.R - event.R_MC)
108  self.hist_theta_residuumhist_theta_residuum.Fill(event.THETA - event.THETA_MC)
109  self.hist_phi_residuumhist_phi_residuum.Fill(event.PHI - event.PHI_MC)
110  self.hist_p_residuumhist_p_residuum.Fill(event.P - event.P_MC)
111 
112  self.hist_chi2hist_chi2.Fill(event.CHI2)
113  if event.R_MC > 1.0:
114  self.hist_chi2_outsidehist_chi2_outside.Fill(event.CHI2)
115  else:
116  assert event.R_MC <= 1.0
117  self.hist_chi2_insidehist_chi2_inside.Fill(event.CHI2)
118 
119  self.hist_mass_vs_mc_masshist_mass_vs_mc_mass.Fill(event.M, event.M_MC)
120 
121  input_root_file.Close()
122  return self
123 
124  @staticmethod
125  def efficiency_plot(found, all, title, x_variable, x_unit, description='', check='', contact='', meta_options=''):
126  """Create an efficiency plot.
127 
128  :param found: Histogram with all found entries (i.e. reconstructed).
129  :param all: Histogram with all entries (i.e. MCTruth).
130  :param title: Title of the histogram.
131  :param x_variable: x variable.
132  :param x_unit: x unit.
133  :param description: Description text shown on the validation page.
134  :param check: Check text shown on the validation page.
135  :param contact: Contact text shown on the validation page.
136  :param meta_options: Meta options for the validation page.
137  :return: ROOT.TEfficiency
138  """
139  efficiency = ROOT.TEfficiency(found, all)
140  efficiency.SetName("".join(title.split()))
141  ylabel = 'Efficiency / ({} {})'.format((found.GetXaxis().GetXmax() -
142  found.GetXaxis().GetXmin()) / found.GetNbinsX(), x_unit)
143  efficiency.SetTitle("{};{} / ({});{}".format(title, x_variable, x_unit, ylabel))
144  efficiency.GetListOfFunctions().Add(ROOT.TNamed('Description', description))
145  efficiency.GetListOfFunctions().Add(ROOT.TNamed('Check', check))
146  efficiency.GetListOfFunctions().Add(ROOT.TNamed('Contact', contact))
147  efficiency.GetListOfFunctions().Add(ROOT.TNamed('MetaOptions', meta_options))
148  return efficiency
149 
150  @staticmethod
151  def histogram_plot(hist, title, x_variable, x_unit=None, description='', check='', contact='', meta_options=''):
152  """Create (annotate) an histogram plot.
153 
154  :param hist: TH1F
155  :param title: Title of the histogram.
156  :param x_variable: x variable.
157  :param x_unit: x unit.
158  :param description: Description text shown on the validation page.
159  :param check: Check text shown on the validation page.
160  :param contact: Contact text shown on the validation page.
161  :param meta_options: Meta options for the validation page.
162  :return: modified hist
163  """
164  hist.SetName("".join(title.split()))
165  xlabel = '{} / ({})'.format(x_variable, x_unit) if x_unit is not None else '{}'.format(x_variable)
166  ylabel = 'Entries / ({} {})'.format((hist.GetXaxis().GetXmax() -
167  hist.GetXaxis().GetXmin()) /
168  hist.GetNbinsX(), x_unit) if x_unit is not None \
169  else 'Entries / ({})'.format((hist.GetXaxis().GetXmax() - hist.GetXaxis().GetXmin()) / hist.GetNbinsX())
170  hist.SetTitle("{};{};{}".format(title, xlabel, ylabel))
171  hist.GetListOfFunctions().Add(ROOT.TNamed('Description', description))
172  hist.GetListOfFunctions().Add(ROOT.TNamed('Check', check))
173  hist.GetListOfFunctions().Add(ROOT.TNamed('Contact', contact))
174  hist.GetListOfFunctions().Add(ROOT.TNamed('MetaOptions', meta_options))
175  return hist
176 
177  @staticmethod
178  def histogram_2d_plot(hist, title, x_variable, y_variable, x_unit=None, y_unit=None,
179  description='', check='', contact='', meta_options=''):
180  """Create a 2d hisogram plot.
181 
182  :param hist: TH2F
183  :param title: Title of the histogram.
184  :param x_variable: x variable
185  :param y_variable: y variable
186  :param x_unit: x unit
187  :param y_unit: y unit
188  :param description: Description text shown on the validation page.
189  :param check: Check text shown on the validation page.
190  :param contact: Contact text shown on the validation page.
191  :param meta_options: Meta options for the validation page.
192  :return: ROOT.TEfficiency
193  :return:
194  """
195  hist.SetName("".join(title.split()))
196  xlabel = '{} / ({})'.format(x_variable, x_unit) if x_unit is not None else '{}'.format(x_variable)
197  ylabel = '{} / ({})'.format(y_variable, y_unit) if y_unit is not None else '{}'.format(y_variable)
198  hist.SetTitle("{};{};{}".format(title, xlabel, ylabel))
199  hist.GetListOfFunctions().Add(ROOT.TNamed('Description', description))
200  hist.GetListOfFunctions().Add(ROOT.TNamed('Check', check))
201  hist.GetListOfFunctions().Add(ROOT.TNamed('Contact', contact))
202  hist.GetListOfFunctions().Add(ROOT.TNamed('MetaOptions', meta_options))
203  return hist
204 
205  def plot(self):
206  """Create plots with the data filled with 'collect_histograms'.
207 
208  :return: self
209  """
210  output_root_file = ROOT.TFile.Open(self.output_fileoutput_file, "RECREATE")
211 
212  V0ValidationPlots.efficiency_plot(self.hist_r_foundhist_r_found, self.hist_rhist_r, 'Efficiency vs R', 'r', 'cm',
213  description='Reconstruction Efficiency vs. r (perpendicular)',
214  check='',
215  contact='software-tracking@belle2.org',
216  meta_options='shifter').Write()
217 
218  V0ValidationPlots.efficiency_plot(self.hist_theta_foundhist_theta_found, self.hist_thetahist_theta, 'Efficiency vs Theta', 'Theta', 'deg',
219  description='Reconstruction Efficiency vs. theta',
220  check='',
221  contact='software-tracking@belle2.org',
222  meta_options='').Write()
223 
224  V0ValidationPlots.efficiency_plot(self.hist_phi_foundhist_phi_found, self.hist_phihist_phi, 'Efficiency vs Phi', 'Phi', 'deg',
225  description='Reconstruction Efficiency vs phi',
226  check='',
227  contact='software-tracking@belle2.org',
228  meta_options='').Write()
229 
230  V0ValidationPlots.efficiency_plot(self.hist_p_foundhist_p_found, self.hist_phist_p, 'Efficiency vs P', 'P', 'GeV',
231  description='Reconstruction Efficiency vs momentum',
232  check='',
233  contact='software-tracking@belle2.org',
234  meta_options='').Write()
235 
236  V0ValidationPlots.histogram_plot(self.hist_invariant_masshist_invariant_mass, "KShort Invariant Mass", "m", "GeV",
237  description='Reconstructed invariant mass of KShorts with standard reconstruction',
238  check='Invariant mass peak around KShort nominal mass 497.61 MeV.',
239  contact='software-tracking@belle2.org',
240  meta_options='shifter').Write()
241 
242  V0ValidationPlots.histogram_plot(self.hist_invariant_mass_residuumhist_invariant_mass_residuum, "KShort Invariant Mass Residuum", "Rec - MC", "GeV",
243  description='Invariant mass residuum of KShorts with standard reconstruction',
244  check='',
245  contact='software-tracking@belle2.org',
246  meta_options='').Write()
247 
248  V0ValidationPlots.histogram_plot(self.hist_chi2hist_chi2, "Chi2 of Vertex Fits.", "Chi2", None,
249  description='Chi2 distributions of the vertex fits.',
250  check='Check if distribution looks like a Chi2 distribution with 1 dof',
251  contact='software-tracking@belle2.org',
252  meta_options='').Write()
253 
254  V0ValidationPlots.histogram_plot(self.hist_chi2_insidehist_chi2_inside, "Chi2 of Vertex Fits Inside Beampipe.", "Chi2", None,
255  description='Chi2 distributions of the vertex fits inside the beampipe.',
256  check='Check if distribution looks like a Chi2 distribution with 1 dof',
257  contact='software-tracking@belle2.org',
258  meta_options='').Write()
259 
260  V0ValidationPlots.histogram_plot(self.hist_chi2_outsidehist_chi2_outside, "Chi2 of Vertex Fits Outside Beampipe.", "Chi2", None,
261  description='Chi2 distributions of the vertex fits outside the beampipe.',
262  check='Check if distribution looks like a Chi2 distribution with 1 dof',
263  contact='software-tracking@belle2.org',
264  meta_options='').Write()
265 
266  V0ValidationPlots.histogram_2d_plot(self.hist_mass_vs_mc_masshist_mass_vs_mc_mass, "Reconstructed vs MC Mass.",
267  "Reconstructed Mass", "GeV", "MC Mass", "GeV",
268  description="Reconstructed mass vs invariant Mass.",
269  check="",
270  contact="software-tracking@belle2.org",
271  meta_options='').Write()
272 
273  V0ValidationPlots.histogram_plot(self.hist_r_residuumhist_r_residuum, "KShort R Residuum", "Rec - MC", "cm",
274  description='R residuum of KShorts with standard reconstruction',
275  check='',
276  contact='software-tracking@belle2.org',
277  meta_options='').Write()
278  V0ValidationPlots.histogram_plot(self.hist_theta_residuumhist_theta_residuum, "KShort Theta Residuum", "Rec - MC", "rad",
279  description='Theta residuum of KShorts with standard reconstruction',
280  check='',
281  contact='software-tracking@belle2.org',
282  meta_options='').Write()
283  V0ValidationPlots.histogram_plot(self.hist_phi_residuumhist_phi_residuum, "KShort Phi Residuum", "Rec - MC", "rad",
284  description='Phi residuum of KShorts with standard reconstruction',
285  check='',
286  contact='software-tracking@belle2.org',
287  meta_options='').Write()
288  V0ValidationPlots.histogram_plot(self.hist_p_residuumhist_p_residuum, "KShort Momentum Residuum", "Rec - MC", "GeV",
289  description='Momentum residuum of KShorts with standard reconstruction',
290  check='',
291  contact='software-tracking@belle2.org',
292  meta_options='').Write()
293 
294  output_root_file.Write()
295  output_root_file.Close()
296  return self
297 
298 
299 if __name__ == '__main__':
hist_chi2_outside
Chi2 of vertex fit outside beampipe.
hist_invariant_mass_res
Invariant mass residual histogram.
def efficiency_plot(found, all, title, x_variable, x_unit, description='', check='', contact='', meta_options='')
hist_mass_vs_mc_mass
2D histogram; invariant mass vs reconstructed mass
def histogram_plot(hist, title, x_variable, x_unit=None, description='', check='', contact='', meta_options='')
hist_chi2_inside
Chi2 of vertex fit inside beampipe.
def histogram_2d_plot(hist, title, x_variable, y_variable, x_unit=None, y_unit=None, description='', check='', contact='', meta_options='')
hist_invariant_mass_residuum
Invariant mass residuum histogram.
def __init__(self, input_file='../V0ValidationHarvested.root', output_file='V0Validation.root')
Definition: plot.py:1