Belle II Software development
2_analysis_pi0.py
1#!/usr/bin/env python3
2
3
10
11"""
12<header>
13 <input>../GenericB_GENSIMRECtoDST.dst.root</input>
14 <output>Pi0_Validation.root</output>
15 <contact>Mario Merola (mario.merola@na.infn.it)</contact>
16 <description>
17 Check the calibration of the ECL in the MC by determining the measured pi0 invariant mass.
18 </description>
19</header>
20"""
21
22import basf2
23import ROOT
24from modularAnalysis import cutAndCopyList, inputMdst
25from stdPi0s import stdPi0s
26from validation_tools.metadata import create_validation_histograms
27from validation_tools.metadata import validation_metadata_update
28from variables import variables as vm
29
30INPUT_FILENAME = "../GenericB_GENSIMRECtoDST.dst.root"
31OUTPUT_FILENAME = "Pi0_Validation.root"
32
33main = basf2.Path()
34inputMdst(INPUT_FILENAME, path=main)
35
36stdPi0s('all', path=main)
37cutAndCopyList('pi0:rec', 'pi0:all', 'daughter(0, E)>0.05 and daughter(1, E)>0.05', path=main)
38cutAndCopyList('pi0:mc', 'pi0:all', 'mcErrors<1', path=main)
39
40vm.addAlias('Mreco', 'M')
41
42
43create_validation_histograms(
44 main, OUTPUT_FILENAME, "pi0:rec",
45 variables_1d=[
46 (
47 "Mreco", 40, 0.08, 0.18,
48 "#pi^{0} reconstructed candidates, invariant mass",
49 "Eldar Ganiev <eldar.ganiev@desy.de>",
50 r"The $\pi^0$ invariant mass distribution with $E_{\gamma}>0.05\, \text{GeV}$",
51 r"Distribution should be peaking at the nominal $\pi^0$ mass.",
52 "M(#pi^{0}) [GeV/c^{2}]", "Candidates", "shifter"
53 ),
54 ],
55 description=r"$\pi^0$ reconstructed mass distribution",
56)
57
58
59vm.addAlias('Mmc', 'M')
60
61create_validation_histograms(
62 main, OUTPUT_FILENAME, "pi0:mc",
63 variables_1d=[
64 (
65 "Mmc", 40, 0.08, 0.18,
66 "#pi^{0} MC candidates, invariant mass",
67 "Eldar Ganiev <eldar.ganiev@desy.de>",
68 r"The $\pi^0$ invariant mass distribution for truth matched candidates",
69 r"Distribution should be peaking at the nominal $\pi^0$ mass.",
70 "M(#pi^{0}) [GeV/c^{2}]", "Candidates", "shifter"
71 ),
72 ],
73 description=r"$\pi^0$ MC mass distribution",
74)
75
76main.add_module('Progress')
77basf2.process(main)
78print(basf2.statistics)
79
80
81f = ROOT.TFile(OUTPUT_FILENAME)
82Mrecohist = f.Get('Mreco')
83Mmchist = f.Get('Mmc')
84
85
86mass = ROOT.RooRealVar("recomass", "m_{#gamma#gamma} [GeV/c^{2}]", 0.11, 0.15)
87
88h_pi0_reco = ROOT.RooDataHist("h_pi0_reco", "h_pi0_reco", ROOT.RooArgList(mass), Mrecohist)
89h_pi0_mc = ROOT.RooDataHist("h_pi0_mc", "h_pi0_mc", ROOT.RooArgList(mass), Mmchist)
90
91
92# pi0 signal PDF is a Gaussian (Crystal Ball also listed in case we want to switch)
93mean = ROOT.RooRealVar("mean", "mean", 0.125, 0.11, 0.15)
94sig1 = ROOT.RooRealVar("#sigma", "sig", 0.007, 0.002, 0.1)
95gau1 = ROOT.RooGaussian("gau1", "gau1", mass, mean, sig1)
96
97alphacb = ROOT.RooRealVar("alphacb", "alpha", 1.5, 0.1, 1.9)
98ncb = ROOT.RooRealVar("ncb", "n", 8) # ,2.,15)
99sigcb = ROOT.RooCBShape("sigcb", "sigcb", mass, mean, sig1, alphacb, ncb)
100
101# pi0 background PDF is a 2nd order Chebyshev
102b1 = ROOT.RooRealVar("b1", "b1", 0.1, -1, 1)
103a1 = ROOT.RooRealVar("a1", "a1", 0.1, -1, 1)
104bList = ROOT.RooArgList(a1, b1)
105bkg = ROOT.RooChebychev("bkg", "bkg", mass, bList)
106
107
108nsig = ROOT.RooRealVar("nsig", "nsig", 3000, 0, 1000000)
109nbkg = ROOT.RooRealVar("nbkg", "nbkg", 12000, 0, 1000000)
110
111
112totalPdf = ROOT.RooAddPdf("totalpdf", "", ROOT.RooArgList(gau1, bkg), ROOT.RooArgList(nsig, nbkg))
113
114
115output = ROOT.TFile("Pi0_Validation_ntuple.root", "recreate")
116
117# Store pi0 mass fit results to a tuple for comparison of mean and width among releases.
118outputNtuple = ROOT.TNtuple(
119 "pi0_mass",
120 "Pi0 mass fit results",
121 "mean:meanerror:width:widtherror:mean_MC:meanerror_MC:width_MC:widtherror_MC")
122
123
124ROOT.gROOT.SetBatch(True)
125canvas = ROOT.TCanvas("canvas", "pi0 mass fit", 1000, 600)
126canvas.Divide(2, 1)
127canvas.cd(1)
128
129# Fit to the reco mass
130totalPdf.fitTo(h_pi0_reco, ROOT.RooFit.Extended(True), ROOT.RooFit.Minimizer("Minuit2", "Migrad"))
131frame1 = mass.frame()
132h_pi0_reco.plotOn(frame1, ROOT.RooFit.Name("Hist"))
133frame1.SetMaximum(frame1.GetMaximum())
134totalPdf.plotOn(frame1, ROOT.RooFit.Name("curve"))
135totalPdf.plotOn(frame1, ROOT.RooFit.Components("gau1"), ROOT.RooFit.LineStyle(ROOT.kDashed), ROOT.RooFit.LineColor(ROOT.kRed))
136totalPdf.plotOn(frame1, ROOT.RooFit.Components("bkg"), ROOT.RooFit.LineStyle(3), ROOT.RooFit.LineColor(ROOT.kBlue))
137frame1.SetMaximum(Mrecohist.GetMaximum() * 1.5)
138frame1.GetXaxis().SetTitleOffset(1.4)
139frame1.GetYaxis().SetTitleOffset(1.5)
140meanval = mean.getVal()
141meanerror = mean.getError()
142width = sig1.getVal()
143widtherror = sig1.getError()
144frame1.Draw("")
145
146canvas.cd(2)
147
148# ---------------------------
149# Fit to the truth matched mass.
150# The same signal parametrisation as for the reco mass but only a 1st order Chebyshev polynomial are used.
151# Re-initialize the fit parameters.
152mean.setVal(0.125)
153sig1.setVal(0.007)
154nsig.setVal(2200)
155nbkg.setVal(700)
156a1.setVal(-0.5)
157
158bkg = ROOT.RooChebychev("bkg1", "bkg", mass, a1)
159totalPdf = ROOT.RooAddPdf("totalpdfMC", "", ROOT.RooArgList(gau1, bkg), ROOT.RooArgList(nsig, nbkg))
160
161# Fit to the truth matched mass
162totalPdf.fitTo(h_pi0_mc, ROOT.RooFit.Extended(True), ROOT.RooFit.Minimizer("Minuit2", "Migrad"))
163frame2 = mass.frame()
164h_pi0_mc.plotOn(frame2, ROOT.RooFit.Name("Hist"))
165frame2.SetMaximum(frame2.GetMaximum())
166totalPdf.plotOn(frame2, ROOT.RooFit.Name("curve"))
167totalPdf.plotOn(frame2, ROOT.RooFit.Components("gau1"), ROOT.RooFit.LineStyle(ROOT.kDashed), ROOT.RooFit.LineColor(ROOT.kRed))
168totalPdf.plotOn(frame2, ROOT.RooFit.Components("bkg1"), ROOT.RooFit.LineStyle(3), ROOT.RooFit.LineColor(ROOT.kBlue))
169frame2.SetMaximum(Mmchist.GetMaximum() * 1.5)
170frame2.GetXaxis().SetTitleOffset(1.4)
171frame2.GetYaxis().SetTitleOffset(1.5)
172meanval_mc = mean.getVal()
173meanerror_mc = mean.getError()
174width_mc = sig1.getVal()
175widtherror_mc = sig1.getError()
176frame2.Draw("")
177
178canvas.Write()
179
180outputNtuple.Fill(meanval, meanerror, width, widtherror, meanval_mc, meanerror_mc, width_mc, widtherror_mc)
181
182outputNtuple.Write()
183
184validation_metadata_update(
185 output,
186 "pi0_mass",
187 title="Pi0 mass fit results",
188 contact="eldar.ganiev@desy.de",
189 description="Fit to the invariant mass of the reconstructed and truth matched pi0s",
190 check="Consistent numerical fit results. Stable mean and width.",
191 metaoptions="shifter")
192
193
194output.Close()