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