23t = ROOT.RooRealVar(
"t",
"time", -1., 15.)
24cosa = ROOT.RooRealVar(
"cosa",
"cos(alpha)", -1., 1.)
28tau = ROOT.RooRealVar(
"tau",
"#tau", 1.5)
29deltaGamma = ROOT.RooRealVar(
"deltaGamma",
"deltaGamma", 0.3)
30tm = ROOT.RooTruthModel(
"tm",
"tm", t)
31coshGBasis = ROOT.RooFormulaVar(
32 "coshGBasis",
"exp(-@0/ @1)*cosh(@0*@2/2)", ROOT.RooArgList(t, tau, deltaGamma))
33sinhGBasis = ROOT.RooFormulaVar(
34 "sinhGBasis",
"exp(-@0/ @1)*sinh(@0*@2/2)", ROOT.RooArgList(t, tau, deltaGamma))
35coshGConv = tm.convolution(coshGBasis, t)
36sinhGConv = tm.convolution(sinhGBasis, t)
39poly1 = ROOT.RooPolyVar(
"poly1",
"poly1", cosa, ROOT.RooArgList(
40 ROOT.RooFit.RooConst(0.5), ROOT.RooFit.RooConst(0.2), ROOT.RooFit.RooConst(0.2)), 0)
41poly2 = ROOT.RooPolyVar(
"poly2",
"poly2", cosa, ROOT.RooArgList(
42 ROOT.RooFit.RooConst(1), ROOT.RooFit.RooConst(-0.2), ROOT.RooFit.RooConst(3)), 0)
45ampl1 = ROOT.RooProduct(
"ampl1",
"amplitude 1",
46 ROOT.RooArgList(poly1, coshGConv))
47ampl2 = ROOT.RooProduct(
"ampl2",
"amplitude 2",
48 ROOT.RooArgList(poly2, sinhGConv))
54f1 = ROOT.RooRealVar(
"f1",
"f1", 1, 0, 2)
55f2 = ROOT.RooRealVar(
"f2",
"f2", 0.5, 0, 2)
58pdf = ROOT.RooRealSumPdf(
"pdf",
"pdf", ROOT.RooArgList(
59 ampl1, ampl2), ROOT.RooArgList(f1, f2))
62data = pdf.generate(ROOT.RooArgSet(t, cosa), 10000)
71hh_cos = ampl1.createHistogram(
"hh_cos", t, ROOT.RooFit.Binning(
72 50), ROOT.RooFit.YVar(cosa, ROOT.RooFit.Binning(50)))
73hh_sin = ampl2.createHistogram(
"hh_sin", t, ROOT.RooFit.Binning(
74 50), ROOT.RooFit.YVar(cosa, ROOT.RooFit.Binning(50)))
75hh_cos.SetLineColor(ROOT.kBlue)
76hh_sin.SetLineColor(ROOT.kRed)
85ras_ampl1 = ROOT.RooArgSet(ampl1)
86pdf.plotOn(frame1, ROOT.RooFit.Components(ras_ampl1),
87 ROOT.RooFit.LineStyle(ROOT.kDashed))
88ras_ampl2 = ROOT.RooArgSet(ampl2)
89pdf.plotOn(frame1, ROOT.RooFit.Components(ras_ampl2), ROOT.RooFit.LineStyle(
90 ROOT.kDashed), ROOT.RooFit.LineColor(ROOT.kRed))
98pdf.plotOn(frame2, ROOT.RooFit.Components(ras_ampl1),
99 ROOT.RooFit.LineStyle(ROOT.kDashed))
100pdf.plotOn(frame2, ROOT.RooFit.Components(ras_ampl2), ROOT.RooFit.LineStyle(
101 ROOT.kDashed), ROOT.RooFit.LineColor(ROOT.kRed))
103c = ROOT.TCanvas(
"rf704_amplitudefit",
"rf704_amplitudefit", 800, 800)
106ROOT.gPad.SetLeftMargin(0.15)
107frame1.GetYaxis().SetTitleOffset(1.4)
110ROOT.gPad.SetLeftMargin(0.15)
111frame2.GetYaxis().SetTitleOffset(1.4)
114ROOT.gPad.SetLeftMargin(0.20)
115hh_cos.GetZaxis().SetTitleOffset(2.3)
118ROOT.gPad.SetLeftMargin(0.20)
119hh_sin.GetZaxis().SetTitleOffset(2.3)
122c.SaveAs(
"rf704_amplitudefit.png")