Belle II Software development
SVDdEdxCalibrationAlgorithm.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#include <svd/calibration/SVDdEdxCalibrationAlgorithm.h>
10#include <svd/dbobjects/SVDdEdxPDFs.h>
11
12#include <TROOT.h>
13#include <TStyle.h>
14#include <TMath.h>
15#include <TFile.h>
16#include <TTree.h>
17#include <TColor.h>
18#include <TLegend.h>
19#include <TCanvas.h>
20#include <TH1D.h>
21#include <TAxis.h>
22
23#include <RooDataSet.h>
24#include <RooRealVar.h>
25#include <RooAddPdf.h>
26#include <RooGaussian.h>
27#include <RooChebychev.h>
28#include <RooBifurGauss.h>
29#include <RooDstD0BG.h>
30#include <RooAbsDataStore.h>
31#include <RooTreeDataStore.h>
32#include <RooMsgService.h>
33#include <RooStats/SPlot.h>
34
35using namespace RooFit;
36using namespace Belle2;
37
39 m_isMakePlots(true)
40{
41 setDescription("SVD dE/dx calibration algorithm");
42}
43
44/* Main calibration method */
46{
47 gROOT->SetBatch(true);
48
49 const auto exprun = getRunList()[0];
50 B2INFO("ExpRun used for calibration: " << exprun.first << " " << exprun.second);
51
52 auto payload = new Belle2::SVDdEdxPDFs();
53
54 // Get data objects
55 auto ttreeLambda = getObjectPtr<TTree>("Lambda");
56 auto ttreeDstar = getObjectPtr<TTree>("Dstar");
57 auto ttreeGamma = getObjectPtr<TTree>("Gamma");
58
59 if (ttreeLambda->GetEntries() < m_MinEvtsPerTree) {
60 B2WARNING("Not enough data for calibration.");
61 return c_NotEnoughData;
62 }
63
64 // call the calibration functions
65 TH2F hLambdaP = LambdaMassFit(ttreeLambda);
66 auto [hDstarK, hDstarPi, hDstarMu] = DstarMassFit(ttreeDstar);
67 TH2F hGammaE = GammaHistogram(ttreeGamma);
68 std::vector<double> pbins = CreatePBinningScheme();
69 TH2F hEmpty("hEmpty", "A histogram returned if we cannot calibrate", m_numPBins, pbins.data(), m_numDEdxBins, 0, m_dedxCutoff);
70 for (int pbin = 0; pbin <= m_numPBins + 1; pbin++) {
71 for (int dedxbin = 0; dedxbin <= m_numDEdxBins + 1; dedxbin++) {
72 hEmpty.SetBinContent(pbin, dedxbin, 0.01);
73 };
74 }
75
76 B2INFO("Histograms are ready, proceed to creating the payload object...");
77 std::vector<TH2F*> hDedxPDFs(6);
78
79 std::array<std::string, 6> part = {"Electron", "Muon", "Pion", "Kaon", "Proton", "Deuteron"};
80
81 TCanvas* candEdx = new TCanvas("candEdx", "SVD dEdx payloads", 1200, 700);
82 candEdx->Divide(3, 2);
83 gStyle->SetOptStat(11);
84
85 for (bool trunmean : {false, true}) {
86 for (int iPart = 0; iPart < 6; iPart++) {
87
88 if (iPart == 0 && trunmean) {
89 hDedxPDFs[iPart] = &hGammaE;
90 payload->setPDF(*hDedxPDFs[iPart], iPart, trunmean);
91 } else if (iPart == 1 && trunmean) {
92 hDedxPDFs[iPart] = &hDstarMu;
93 payload->setPDF(*hDedxPDFs[iPart], iPart, trunmean);
94 } else if (iPart == 2 && trunmean) {
95 hDedxPDFs[iPart] = &hDstarPi;
96 payload->setPDF(*hDedxPDFs[iPart], iPart, trunmean);
97 } else if (iPart == 3 && trunmean) {
98 hDedxPDFs[iPart] = &hDstarK;
99 payload->setPDF(*hDedxPDFs[iPart], iPart, trunmean);
100 } else if (iPart == 4 && trunmean) {
101 hDedxPDFs[iPart] = &hLambdaP;
102 payload->setPDF(*hDedxPDFs[iPart], iPart, trunmean);
103 } else if (iPart == 5 && trunmean) {
104 hDedxPDFs[iPart] = &hEmpty;
105 hDedxPDFs[iPart]->SetName("hist_d1_1000010020_trunc");
106 payload->setPDF(*hDedxPDFs[iPart], iPart, trunmean);
107 } else if (iPart == 0 && !trunmean) {
108 hDedxPDFs[iPart] = &hEmpty;
109 hDedxPDFs[iPart]->SetName("hist_d1_11");
110 payload->setPDF(*hDedxPDFs[iPart], iPart, trunmean);
111 } else if (iPart == 1 && !trunmean) {
112 hDedxPDFs[iPart] = &hEmpty;
113 hDedxPDFs[iPart]->SetName("hist_d1_13");
114 payload->setPDF(*hDedxPDFs[iPart], iPart, trunmean);
115 } else if (iPart == 2 && !trunmean) {
116 hDedxPDFs[iPart] = &hEmpty;
117 hDedxPDFs[iPart]->SetName("hist_d1_211");
118 payload->setPDF(*hDedxPDFs[iPart], iPart, trunmean);
119 } else if (iPart == 3 && !trunmean) {
120 hDedxPDFs[iPart] = &hEmpty;
121 hDedxPDFs[iPart]->SetName("hist_d1_321");
122 payload->setPDF(*hDedxPDFs[iPart], iPart, trunmean);
123 } else if (iPart == 4 && !trunmean) {
124 hDedxPDFs[iPart] = &hEmpty;
125 hDedxPDFs[iPart]->SetName("hist_d1_2212");
126 payload->setPDF(*hDedxPDFs[iPart], iPart, trunmean);
127 } else if (iPart == 5 && !trunmean) {
128 hDedxPDFs[iPart] = &hEmpty;
129 hDedxPDFs[iPart]->SetName("hist_d1_1000010020");
130 payload->setPDF(*hDedxPDFs[iPart], iPart, trunmean);
131 }
132
133 else
134 hDedxPDFs[iPart] = &hEmpty;
135 payload->setPDF(*hDedxPDFs[iPart], iPart, trunmean);
136
137 candEdx->cd(iPart + 1);
138 hDedxPDFs[iPart]->SetTitle(Form("%s; p(GeV/c) of %s; dE/dx", hDedxPDFs[iPart]->GetTitle(), part[iPart].data()));
139 hDedxPDFs[iPart]->DrawCopy("colz");
140
141 if (m_isMakePlots) {
142 candEdx->SaveAs("PlotsSVDdEdxPDFs_wTruncMean.pdf");
143 }
144 }
145 // candEdx->SetTitle(Form("Likehood dist. of charged particles from %s, trunmean = %s", idet.data(), check.str().data()));
146 }
147
148 saveCalibration(payload, "SVDdEdxPDFs");
149 B2INFO("SVD dE/dx calibration done!");
150
151 return c_OK;
152}
153
154TH2F SVDdEdxCalibrationAlgorithm::LambdaMassFit(std::shared_ptr<TTree> preselTree)
155{
156 B2INFO("Configuring the Lambda fit...");
157 gROOT->SetBatch(true);
158 RooMsgService::instance().setGlobalKillBelow(RooFit::WARNING);
159
160 RooRealVar InvM("InvM", "m(p^{+}#pi^{-})", 1.1, 1.13, "GeV/c^{2}");
161
162 RooRealVar ProtonMomentum("ProtonMomentum", "momentum for p", -1.e8, 1.e8);
163 RooRealVar ProtonSVDdEdx("ProtonSVDdEdx", "", -1.e8, 1.e8);
164
165 RooRealVar exp("exp", "experiment number", 0, 1.e5);
166 RooRealVar run("run", "run number", 0, 1.e7);
167
168 auto variables = new RooArgSet();
169
170 variables->add(InvM);
171
172 variables->add(ProtonMomentum);
173 variables->add(ProtonSVDdEdx);
174 variables->add(exp);
175 variables->add(run);
176
177 RooDataSet* LambdaDataset = new RooDataSet("LambdaDataset", "LambdaDataset", preselTree.get(), *variables);
178
179 if (LambdaDataset->sumEntries() == 0) {
180 B2FATAL("The Lambda dataset is empty, stopping here");
181 }
182
183 // the signal PDF; might be revisited at a later point
184
185 RooRealVar GaussMean("GaussMean", " GaussMean", 1.116, 1.111, 1.12);
186 RooRealVar GaussSigma("GaussSigma", "#sigma_{1}", 3.e-3, 3.e-5, 10.e-3);
187 RooGaussian LambdaGauss("LambdaGauss", "LambdaGauss", InvM, GaussMean, GaussSigma);
188
189 /* temporary RooRealVar sigmaBifurGaussL1 and sigmaBifurGaussR1 to replace
190 * RooRealVar resolutionParamL("resolutionParamL", "resolutionParamL", 0.4, 5.e-4, 1.0);
191 * RooRealVar resolutionParamR("resolutionParamR", "resolutionParamR", 0.4, 5.e-4, 1.0);
192 * RooFormulaVar sigmaBifurGaussL1("sigmaBifurGaussL1", "resolutionParamL*GaussSigma", RooArgSet(resolutionParamL, GaussSigma));
193 * RooFormulaVar sigmaBifurGaussR1("sigmaBifurGaussR1", "resolutionParamR*GaussSigma", RooArgSet(resolutionParamR, GaussSigma));
194 */
195 RooRealVar sigmaBifurGaussL1("sigmaBifurGaussL1", "sigma left", 0.4 * 3.e-3, 3.e-5, 10.e-3);
196 RooRealVar sigmaBifurGaussR1("sigmaBifurGaussR1", "sigma right", 0.4 * 3.e-3, 3.e-5, 10.e-3);
197 RooBifurGauss LambdaBifurGauss("LambdaBifurGauss", "LambdaBifurGauss", InvM, GaussMean, sigmaBifurGaussL1, sigmaBifurGaussR1);
198
199 /* temporary RooRealVar sigmaBifurGaussL2 to replace
200 * RooRealVar resolutionParam2("resolutionParam2", "resolutionParam2", 0.2, 5.e-4, 1.0);
201 * sigmaBifurGaussL2("sigmaBifurGaussL2", "resolutionParam2*GaussSigma", RooArgSet(resolutionParam2, GaussSigma));
202 */
203 RooRealVar sigmaBifurGaussL2("sigmaBifurGaussL2", "sigmaBifurGaussL2", 0.2 * 3.e-3, 3.e-5, 10.e-3);
204 RooGaussian LambdaBifurGauss2("LambdaBifurGauss2", "LambdaBifurGauss2", InvM, GaussMean, sigmaBifurGaussL2);
205
206 RooRealVar fracBifurGaussYield("fracBifurGaussYield", "fracBifurGaussYield", 0.3, 5.e-4, 1.0);
207 RooRealVar fracGaussYield("fracGaussYield", "fracGaussYield", 0.8, 5.e-4, 1.0);
208
209 RooAddPdf LambdaCombinedBifurGauss("LambdaCombinedBifurGauss", "LambdaBifurGauss + LambdaBifurGauss2 ", RooArgList(LambdaBifurGauss,
210 LambdaBifurGauss2), RooArgList(fracBifurGaussYield));
211
212 RooAddPdf LambdaSignalPDF("LambdaSignalPDF", "LambdaCombinedBifurGauss + LambdaGauss", RooArgList(LambdaCombinedBifurGauss,
213 LambdaGauss), RooArgList(fracGaussYield));
214
215 // Background PDF
216 RooRealVar BkgPolyCoef0("BkgPolyCoef0", "BkgPolyCoef0", 0.1, 0., 1.5);
217 RooRealVar BkgPolyCoef1("BkgPolyCoef1", "BkgPolyCoef1", -0.5, -1.5, -1.e-3);
218 RooChebychev BkgPolyPDF("BkgPolyPDF", "BkgPolyPDF", InvM, RooArgList(BkgPolyCoef0, BkgPolyCoef1));
219
220 RooRealVar nSignalLambda("nSignalLambda", "nSignalLambda", 0.6 * preselTree->GetEntries(), 0., 0.99 * preselTree->GetEntries());
221 RooRealVar nBkgLambda("nBkgLambda", "nBkgLambda", 0.4 * preselTree->GetEntries(), 0., 0.99 * preselTree->GetEntries());
222 RooAddPdf totalPDFLambda("totalPDFLambda", "totalPDFLambda pdf", RooArgList(LambdaSignalPDF, BkgPolyPDF),
223 RooArgList(nSignalLambda, nBkgLambda));
224
225 B2INFO("Lambda: Start fitting...");
226 RooFitResult* LambdaFitResult = totalPDFLambda.fitTo(*LambdaDataset, Save(kTRUE), PrintLevel(-1));
227
228 int status = LambdaFitResult->status();
229 int covqual = LambdaFitResult->covQual();
230 double diff = nSignalLambda.getValV() + nBkgLambda.getValV() - LambdaDataset->sumEntries();
231
232 B2INFO("Lambda: Fit status: " << status << "; covariance quality: " << covqual);
233 // if the fit is not healthy, try again once before giving up, with a slightly different setup:
234 if ((status > 0) || (TMath::Abs(diff) > 1.) || (nSignalLambda.getError() < sqrt(nSignalLambda.getValV()))
235 || (nSignalLambda.getError() > (nSignalLambda.getValV()))) {
236
237 LambdaFitResult = totalPDFLambda.fitTo(*LambdaDataset, Save(), Strategy(2), Offset(1));
238 status = LambdaFitResult->status();
239 covqual = LambdaFitResult->covQual();
240 diff = nSignalLambda.getValV() + nBkgLambda.getValV() - LambdaDataset->sumEntries();
241 }
242
243 if ((status > 0) || (TMath::Abs(diff) > 1.) || (nSignalLambda.getError() < sqrt(nSignalLambda.getValV()))
244 || (nSignalLambda.getError() > (nSignalLambda.getValV()))) {
245 B2WARNING("Lambda: Fit problem: fit status " << status << "; sum of component yields minus the dataset yield is " << diff <<
246 "; signal yield is " << nSignalLambda.getValV() << ", while its uncertainty is " << nSignalLambda.getError());
247 }
248 if (covqual < 2) {
249 B2INFO("Lambda: Fit warning: covariance quality " << covqual);
250 }
251
252 TCanvas* canvLambda = new TCanvas("canvLambda", "canvLambda");
253 RooPlot* LambdaFitFrame = LambdaDataset->plotOn(InvM.frame(130));
254 totalPDFLambda.plotOn(LambdaFitFrame, LineColor(TColor::GetColor("#4575b4")));
255
256 double chisquare = LambdaFitFrame->chiSquare();
257 B2INFO("Lambda: Fit chi2 = " << chisquare);
258 totalPDFLambda.paramOn(LambdaFitFrame, Layout(0.6, 0.96, 0.93), Format("NEU", AutoPrecision(2)));
259 LambdaFitFrame->getAttText()->SetTextSize(0.03);
260
261 totalPDFLambda.plotOn(LambdaFitFrame, Components("LambdaSignalPDF"), LineColor(TColor::GetColor("#d73027")));
262 totalPDFLambda.plotOn(LambdaFitFrame, Components("BkgPolyPDF"), LineColor(TColor::GetColor("#fc8d59")));
263 totalPDFLambda.plotOn(LambdaFitFrame, LineColor(TColor::GetColor("#4575b4")));
264
265 LambdaFitFrame->GetXaxis()->SetTitle("m(p#pi^{-}) (GeV/c^{2})");
266
267 LambdaFitFrame->Draw();
268
269 if (m_isMakePlots) {
270 canvLambda->Print("SVDdEdxCalibrationFitLambda.pdf");
271 TFile LambdaFitPlotFile("SVDdEdxCalibrationLambdaFitPlotFile.root", "RECREATE");
272 canvLambda->Write();
273 LambdaFitPlotFile.Close();
274 }
275 RooStats::SPlot* sPlotDatasetLambda = new RooStats::SPlot("sData", "An SPlot", *LambdaDataset, &totalPDFLambda,
276 RooArgList(nSignalLambda, nBkgLambda));
277
278 for (int iEvt = 0; iEvt < 5; iEvt++) {
279 if (TMath::Abs(sPlotDatasetLambda->GetSWeight(iEvt, "nSignalLambda") + sPlotDatasetLambda->GetSWeight(iEvt,
280 "nBkgLambda") - 1) > 5.e-3)
281 B2FATAL("Lambda: sPlot error: sum of weights not equal to 1");
282 }
283
284 RooDataSet* LambdaDatasetSWeighted = new RooDataSet(LambdaDataset->GetName(), LambdaDataset->GetTitle(), LambdaDataset,
285 *LambdaDataset->get());
286
287 TTree* treeLambda_sw = LambdaDatasetSWeighted->GetClonedTree();
288 treeLambda_sw->SetName("treeLambda_sw");
289
290 B2INFO("Lambda: sPlot done. Proceed to histogramming");
291
292 std::vector<double> pbins = CreatePBinningScheme();
293
294 TH2F* hLambdaP = new TH2F("hist_d1_2212_trunc", "hist_d1_2212_trunc", m_numPBins, pbins.data(), m_numDEdxBins, 0, m_dedxCutoff);
295
296 treeLambda_sw->Draw("ProtonSVDdEdx:ProtonMomentum>>hist_d1_2212_trunc",
297 "nSignalLambda_sw * (ProtonMomentum>0.15) * (ProtonSVDdEdx>0)", "goff");
298
299 // produce the 1D profile (for data-MC comparisons)
300 if (m_isMakePlots) {
301 TH1F* ProtonProfile = (TH1F*)hLambdaP->ProfileX("ProtonProfile");
302 ProtonProfile->SetTitle("ProtonProfile");
303 canvLambda->SetTicky(1);
304 ProtonProfile->GetYaxis()->SetRangeUser(0, m_dedxCutoff);
305 ProtonProfile->GetXaxis()->SetTitle("Momentum, GeV/c");
306 ProtonProfile->GetYaxis()->SetTitle("dE/dx");
307 ProtonProfile->Draw();
308 canvLambda->Print("SVDdEdxCalibrationProfileProton.pdf");
309 TFile ProtonProfileFile("SVDdEdxCalibrationProfileProton.root", "RECREATE");
310 ProtonProfile->Write();
311 ProtonProfileFile.Close();
312 canvLambda->SetTicky(0);
313 }
314
315 // for each momentum bin, normalize the pdf
316
317 // hLambdaP normalisation
318 for (int pbin = 0; pbin <= m_numPBins + 1; pbin++) {
319 for (int dedxbin = 0; dedxbin <= m_numDEdxBins + 1; dedxbin++) {
320 // get rid of the bins with negative weights
321 if (hLambdaP->GetBinContent(pbin, dedxbin) <= 1) {
322 hLambdaP->SetBinContent(pbin, dedxbin, 0);
323 };
324 }
325 // create a projection (1D histogram) in a given momentum bin
326 TH1D* slice = (TH1D*)hLambdaP->ProjectionY("slice", pbin, pbin);
327 // normalise, but ignore the cases with empty histograms
328 if (slice->Integral() > 0) {
329 slice->Scale(1. / slice->Integral());
330 }
331 // fill back the 2D histo with the result
332 for (int dedxbin = 0; dedxbin <= m_numDEdxBins + 1; dedxbin++) {
333 hLambdaP->SetBinContent(pbin, dedxbin, slice->GetBinContent(dedxbin));
334 }
335 }
336
337 if (m_isMakePlots) {
338 hLambdaP->Draw("COLZ");
339 canvLambda->Print("SVDdEdxCalibrationHistoLambda.pdf");
340 }
341
342 return *hLambdaP;
343}
344
345std::tuple<TH2F, TH2F, TH2F> SVDdEdxCalibrationAlgorithm::DstarMassFit(std::shared_ptr<TTree> preselTree)
346{
347 B2INFO("Configuring the Dstar fit...");
348 gROOT->SetBatch(true);
349 RooMsgService::instance().setGlobalKillBelow(RooFit::WARNING);
350
351 RooRealVar deltaM("deltaM", "m(D*)-m(D^{0})", 0.139545, 0.151, "GeV/c^{2}");
352
353 RooRealVar KaonMomentum("KaonMomentum", "momentum for Kaon(GeV)", -1.e8, 1.e8);
354 RooRealVar KaonSVDdEdx("KaonSVDdEdx", "", -1.e8, 1.e8);
355 RooRealVar PionDMomentum("PionDMomentum", "momentum for pion(GeV)", -1.e8, 1.e8);
356 RooRealVar PionDSVDdEdx("PionDSVDdEdx", "", -1.e8, 1.e8);
357 RooRealVar SlowPionMomentum("SlowPionMomentum", "momentum for slow pion(GeV)", -1.e8, 1.e8);
358 RooRealVar SlowPionSVDdEdx("SlowPionSVDdEdx", "", -1.e8, 1.e8);
359
360 RooRealVar exp("exp", "experiment number", 0, 1.e5);
361 RooRealVar run("run", "run number", 0, 1.e8);
362
363 auto variables = new RooArgSet();
364 variables->add(deltaM);
365 variables->add(KaonMomentum);
366 variables->add(KaonSVDdEdx);
367 variables->add(PionDMomentum);
368 variables->add(PionDSVDdEdx);
369 variables->add(SlowPionMomentum);
370 variables->add(SlowPionSVDdEdx);
371 variables->add(exp);
372 variables->add(run);
373
374 RooDataSet* DstarDataset = new RooDataSet("DstarDataset", "DstarDataset", preselTree.get(), *variables);
375
376 if (DstarDataset->sumEntries() == 0) {
377 B2FATAL("The Dstar dataset is empty, stopping here");
378 }
379
380 RooPlot* DstarFitFrame = DstarDataset->plotOn(deltaM.frame());
381
382 RooRealVar GaussMean("GaussMean", "GaussMean", 0.145, 0.140, 0.150);
383 RooRealVar GaussSigma1("GaussSigma1", "GaussSigma1", 0.01, 1.e-4, 1.0);
384 RooGaussian DstarGauss1("DstarGauss1", "DstarGauss1", deltaM, GaussMean, GaussSigma1);
385 RooRealVar GaussSigma2("GaussSigma2", "GaussSigma2", 0.001, 1.e-4, 1.0);
386 RooGaussian DstarGauss2("DstarGauss2", "DstarGauss2", deltaM, GaussMean, GaussSigma2);
387 RooRealVar fracGaussYield("fracGaussYield", "Fraction of two Gaussians", 0.75, 0.0, 1.0);
388 RooAddPdf DstarSignalPDF("DstarSignalPDF", "DstarGauss1+DstarGauss2", RooArgList(DstarGauss1, DstarGauss2), fracGaussYield);
389
390 RooRealVar dm0Bkg("dm0Bkg", "dm0", 0.13957018, 0.130, 0.140);
391 RooRealVar aBkg("aBkg", "a", -0.0784, -0.08, 3.0);
392 RooRealVar bBkg("bBkg", "b", -0.444713, -0.5, 0.4);
393 RooRealVar cBkg("cBkg", "c", 0.3);
394 RooDstD0BG DstarBkgPDF("DstarBkgPDF", "DstarBkgPDF", deltaM, dm0Bkg, cBkg, aBkg, bBkg);
395 RooRealVar nSignalDstar("nSignalDstar", "signal yield", 0.5 * preselTree->GetEntries(), 0, preselTree->GetEntries());
396 RooRealVar nBkgDstar("nBkgDstar", "background yield", 0.5 * preselTree->GetEntries(), 0, preselTree->GetEntries());
397 RooAddPdf totalPDFDstar("totalPDFDstar", "totalPDFDstar pdf", RooArgList(DstarSignalPDF, DstarBkgPDF),
398 RooArgList(nSignalDstar, nBkgDstar));
399
400 B2INFO("Dstar: Start fitting...");
401 RooFitResult* DstarFitResult = totalPDFDstar.fitTo(*DstarDataset, Save(kTRUE), PrintLevel(-1));
402
403 int status = DstarFitResult->status();
404 int covqual = DstarFitResult->covQual();
405 double diff = nSignalDstar.getValV() + nBkgDstar.getValV() - DstarDataset->sumEntries();
406
407 B2INFO("Dstar: Fit status: " << status << "; covariance quality: " << covqual);
408 // if the fit is not healthy, try again once before giving up, with a slightly different setup:
409 if ((status > 0) || (TMath::Abs(diff) > 1.) || (nSignalDstar.getError() < sqrt(nSignalDstar.getValV()))
410 || (nSignalDstar.getError() > (nSignalDstar.getValV()))) {
411
412 DstarFitResult = totalPDFDstar.fitTo(*DstarDataset, Save(), Strategy(2), Offset(1));
413 status = DstarFitResult->status();
414 covqual = DstarFitResult->covQual();
415 diff = nSignalDstar.getValV() + nBkgDstar.getValV() - DstarDataset->sumEntries();
416 }
417
418 if ((status > 0) || (TMath::Abs(diff) > 1.) || (nSignalDstar.getError() < sqrt(nSignalDstar.getValV()))
419 || (nSignalDstar.getError() > (nSignalDstar.getValV()))) {
420 B2WARNING("Dstar: Fit problem: fit status " << status << "; sum of component yields minus the dataset yield is " << diff <<
421 "; signal yield is " << nSignalDstar.getValV() << ", while its uncertainty is " << nSignalDstar.getError());
422 }
423 if (covqual < 2) {
424 B2INFO("Dstar: Fit warning: covariance quality " << covqual);
425 }
426
427 totalPDFDstar.plotOn(DstarFitFrame, LineColor(TColor::GetColor("#4575b4")));
428
429 double chisquare = DstarFitFrame->chiSquare();
430 B2INFO("Dstar: Fit chi2 = " << chisquare);
431 totalPDFDstar.paramOn(DstarFitFrame, Layout(0.63, 0.96, 0.93), Format("NEU", AutoPrecision(2)));
432 DstarFitFrame->getAttText()->SetTextSize(0.03);
433
434 totalPDFDstar.plotOn(DstarFitFrame, Components("DstarSignalPDF"), LineColor(TColor::GetColor("#d73027")));
435 totalPDFDstar.plotOn(DstarFitFrame, Components("DstarBkgPDF"), LineColor(TColor::GetColor("#fc8d59")));
436 totalPDFDstar.plotOn(DstarFitFrame, LineColor(TColor::GetColor("#4575b4")));
437
438 DstarFitFrame->GetXaxis()->SetTitle("#Deltam [GeV/c^{2}]");
439 TCanvas* canvDstar = new TCanvas("canvDstar", "canvDstar");
440 canvDstar->cd();
441
442 DstarFitFrame->Draw();
443
444 if (m_isMakePlots) {
445 canvDstar->Print("SVDdEdxCalibrationFitDstar.pdf");
446 TFile DstarFitPlotFile("SVDdEdxCalibrationDstarFitPlotFile.root", "RECREATE");
447 canvDstar->Write();
448 DstarFitPlotFile.Close();
449 }
450
452
453 RooStats::SPlot* sPlotDatasetDstar = new RooStats::SPlot("sData", "An SPlot", *DstarDataset, &totalPDFDstar,
454 RooArgList(nSignalDstar, nBkgDstar));
455
456 for (int iEvt = 0; iEvt < 5; iEvt++) {
457 if (TMath::Abs(sPlotDatasetDstar->GetSWeight(iEvt, "nSignalDstar") + sPlotDatasetDstar->GetSWeight(iEvt, "nBkgDstar") - 1) > 5.e-3)
458 B2FATAL("Dstar: sPlot error: sum of weights not equal to 1");
459 }
460
461 RooDataSet* DstarDatasetSWeighted = new RooDataSet(DstarDataset->GetName(), DstarDataset->GetTitle(), DstarDataset,
462 *DstarDataset->get());
463
464 TTree* treeDstar_sw = DstarDatasetSWeighted->GetClonedTree();
465 treeDstar_sw->SetName("treeDstar_sw");
466
467 B2INFO("Dstar: sPlot done. Proceed to histogramming");
468
469 std::vector<double> pbins = CreatePBinningScheme();
470
471 // the kaon payload
472 TH2F* hDstarK = new TH2F("hist_d1_321_trunc", "hist_d1_321_trunc", m_numPBins, pbins.data(),
474 // the pion payload
475 TH2F* hDstarPi = new TH2F("hist_d1_211_trunc", "hist_d1_211_trunc", m_numPBins, pbins.data(),
477
478 treeDstar_sw->Draw("KaonSVDdEdx:KaonMomentum>>hist_d1_321_trunc", "nSignalDstar_sw * (KaonSVDdEdx>0)", "goff");
479 // the pion one will be built from both pions in the Dstar decay tree
480 TH2F* hDstarPiPart1 = (TH2F*)hDstarPi->Clone("hist_d1_211_truncPart1");
481 TH2F* hDstarPiPart2 = (TH2F*)hDstarPi->Clone("hist_d1_211_truncPart2");
482
483 treeDstar_sw->Draw("PionDSVDdEdx:PionDMomentum>>hist_d1_211_truncPart1", "nSignalDstar_sw * (PionDSVDdEdx>0)", "goff");
484 treeDstar_sw->Draw("SlowPionSVDdEdx:SlowPionMomentum>>hist_d1_211_truncPart2", "nSignalDstar_sw * (SlowPionSVDdEdx>0)", "goff");
485 hDstarPi->Add(hDstarPiPart1);
486 hDstarPi->Add(hDstarPiPart2);
487
488 // the current strategy assumes that the muon and pion payloads are indistinguishable: clone the pion one
489 TH2F* hDstarMu = (TH2F*)hDstarPi->Clone("hist_d1_13_trunc");
490 hDstarMu->SetTitle("hist_d1_13_trunc");
491
492 // produce the 1D profile (for data-MC comparisons)
493 if (m_isMakePlots) {
494 TH1F* PionProfile = (TH1F*)hDstarPi->ProfileX("PionProfile");
495 PionProfile->SetTitle("PionProfile");
496 canvDstar->SetTicky(1);
497 PionProfile->GetYaxis()->SetRangeUser(0, m_dedxCutoff);
498 PionProfile->GetXaxis()->SetTitle("Momentum, GeV/c");
499 PionProfile->GetYaxis()->SetTitle("dE/dx");
500 PionProfile->Draw();
501 canvDstar->Print("SVDdEdxCalibrationProfilePion.pdf");
502 TFile PionProfileFile("SVDdEdxCalibrationProfilePion.root", "RECREATE");
503 PionProfile->Write();
504 PionProfileFile.Close();
505
506 TH1F* KaonProfile = (TH1F*)hDstarK->ProfileX("KaonProfile");
507 KaonProfile->SetTitle("KaonProfile");
508 KaonProfile->GetYaxis()->SetRangeUser(0, m_dedxCutoff);
509 KaonProfile->GetXaxis()->SetTitle("Momentum, GeV/c");
510 KaonProfile->GetYaxis()->SetTitle("dE/dx");
511 KaonProfile->Draw();
512 canvDstar->Print("SVDdEdxCalibrationProfileKaon.pdf");
513 TFile KaonProfileFile("SVDdEdxCalibrationProfileKaon.root", "RECREATE");
514 KaonProfile->Write();
515 KaonProfileFile.Close();
516 canvDstar->SetTicky(0);
517 }
518
519 // hDstarK normalisation
520 // for each momentum bin, normalize the pdf
521
522 for (int pbin = 0; pbin <= m_numPBins + 1; pbin++) {
523 for (int dedxbin = 0; dedxbin <= m_numDEdxBins + 1; dedxbin++) {
524 // get rid of the bins with negative weights
525 if (hDstarK->GetBinContent(pbin, dedxbin) <= 1) {
526 hDstarK->SetBinContent(pbin, dedxbin, 0);
527 };
528 }
529 // create a projection (1D histogram) in a given momentum bin
530 TH1D* slice = (TH1D*)hDstarK->ProjectionY("slice", pbin, pbin);
531 // normalise, but ignore the cases with empty histograms
532 if (slice->Integral() > 0) {
533 slice->Scale(1. / slice->Integral());
534 }
535 // fill back the 2D histo with the result
536 for (int dedxbin = 0; dedxbin <= m_numDEdxBins + 1; dedxbin++) {
537 hDstarK->SetBinContent(pbin, dedxbin, slice->GetBinContent(dedxbin));
538 }
539 }
540
541 // hDstarPi normalisation
542 for (int pbin = 0; pbin <= m_numPBins + 1; pbin++) {
543 for (int dedxbin = 0; dedxbin <= m_numDEdxBins + 1; dedxbin++) {
544 // get rid of the bins with negative weights
545 if (hDstarPi->GetBinContent(pbin, dedxbin) <= 1) {
546 hDstarPi->SetBinContent(pbin, dedxbin, 0);
547 };
548 }
549 // create a projection (1D histogram) in a given momentum bin
550 TH1D* slice = (TH1D*)hDstarPi->ProjectionY("slice", pbin, pbin);
551 // normalise, but ignore the cases with empty histograms
552 if (slice->Integral() > 0) {
553 slice->Scale(1. / slice->Integral());
554 }
555 // fill back the 2D histo with the result
556 for (int dedxbin = 0; dedxbin <= m_numDEdxBins + 1; dedxbin++) {
557 hDstarPi->SetBinContent(pbin, dedxbin, slice->GetBinContent(dedxbin));
558 }
559 }
560
561 // hDstarMu normalisation
562 for (int pbin = 0; pbin <= m_numPBins + 1; pbin++) {
563 for (int dedxbin = 0; dedxbin <= m_numDEdxBins + 1; dedxbin++) {
564 // get rid of the bins with negative weights
565 if (hDstarMu->GetBinContent(pbin, dedxbin) <= 1) {
566 hDstarMu->SetBinContent(pbin, dedxbin, 0);
567 };
568 }
569 // create a projection (1D histogram) in a given momentum bin
570 TH1D* slice = (TH1D*)hDstarMu->ProjectionY("slice", pbin, pbin);
571 // normalise, but ignore the cases with empty histograms
572 if (slice->Integral() > 0) {
573 slice->Scale(1. / slice->Integral());
574 }
575 // fill back the 2D histo with the result
576 for (int dedxbin = 0; dedxbin <= m_numDEdxBins + 1; dedxbin++) {
577 hDstarMu->SetBinContent(pbin, dedxbin, slice->GetBinContent(dedxbin));
578 }
579 }
580 if (m_isMakePlots) {
581 hDstarK->Draw("COLZ");
582 canvDstar->Print("SVDdEdxCalibrationHistoDstarK.pdf");
583 hDstarPi->Draw("COLZ");
584 canvDstar->Print("SVDdEdxCalibrationHistoDstarPi.pdf");
585 hDstarMu->Draw("COLZ");
586 canvDstar->Print("SVDdEdxCalibrationHistoDstarMu.pdf");
587 }
588
589 return std::make_tuple(*hDstarK, *hDstarPi, *hDstarMu);
590}
591
592TH2F SVDdEdxCalibrationAlgorithm::GammaHistogram(std::shared_ptr<TTree> preselTree)
593{
594 B2INFO("Histogramming the converted photon selection...");
595 gROOT->SetBatch(true);
596
597 if (preselTree->GetEntries() == 0) {
598 B2FATAL("The Gamma tree is empty, stopping here");
599 }
600 std::vector<double> pbins = CreatePBinningScheme();
601
602 TH2F* hGammaE = new TH2F("hist_d1_11_trunc", "hist_d1_11_trunc", m_numPBins, pbins.data(), m_numDEdxBins, 0, m_dedxCutoff);
603
604 TH2F* hGammaEPart1 = (TH2F*)hGammaE->Clone("hist_d1_11_truncPart1");
605 TH2F* hGammaEPart2 = (TH2F*)hGammaE->Clone("hist_d1_11_truncPart2");
606
607 preselTree->Draw("FirstElectronSVDdEdx:FirstElectronMomentum>>hist_d1_11_truncPart1", "FirstElectronSVDdEdx>0", "goff");
608 preselTree->Draw("SecondElectronSVDdEdx:SecondElectronMomentum>>hist_d1_11_truncPart2", "SecondElectronSVDdEdx>0", "goff");
609 hGammaE->Add(hGammaEPart1);
610 hGammaE->Add(hGammaEPart2);
611
612 // produce the 1D profile (for data-MC comparisons)
613 TCanvas* canvGamma = new TCanvas("canvGamma", "canvGamma");
614 if (m_isMakePlots) {
615 TH1F* ElectronProfile = (TH1F*)hGammaE->ProfileX("ElectronProfile");
616 ElectronProfile->SetTitle("ElectronProfile");
617 canvGamma->SetTicky(1);
618 ElectronProfile->GetYaxis()->SetRangeUser(0, m_dedxCutoff);
619 ElectronProfile->GetXaxis()->SetTitle("Momentum, GeV/c");
620 ElectronProfile->GetYaxis()->SetTitle("dE/dx");
621 ElectronProfile->Draw();
622 canvGamma->Print("SVDdEdxCalibrationProfileElectron.pdf");
623 TFile ElectronProfileFile("SVDdEdxCalibrationProfileElectron.root", "RECREATE");
624 ElectronProfile->Write();
625 ElectronProfileFile.Close();
626 canvGamma->SetTicky(0);
627 }
628
629 // for each momentum bin, normalize the pdf
630 // hGammaE normalisation
631 for (int pbin = 0; pbin <= m_numPBins + 1; pbin++) {
632 for (int dedxbin = 0; dedxbin <= m_numDEdxBins + 1; dedxbin++) {
633 // get rid of the bins with negative weights
634 if (hGammaE->GetBinContent(pbin, dedxbin) <= 1) {
635 hGammaE->SetBinContent(pbin, dedxbin, 0);
636 };
637 }
638
639 // create a projection (1D histogram) in a given momentum bin
640 TH1D* slice = (TH1D*)hGammaE->ProjectionY("slice", pbin, pbin);
641 // normalise, but ignore the cases with empty histograms
642 if (slice->Integral() > 0) {
643 slice->Scale(1. / slice->Integral());
644 }
645 // fill back the 2D histo with the result
646 for (int dedxbin = 0; dedxbin <= m_numDEdxBins + 1; dedxbin++) {
647 hGammaE->SetBinContent(pbin, dedxbin, slice->GetBinContent(dedxbin));
648 }
649 }
650
651 if (m_isMakePlots) {
652 hGammaE->Draw("COLZ");
653 canvGamma->Print("SVDdEdxCalibrationHistoGamma.pdf");
654 }
655
656 return *hGammaE;
657}
Base class for calibration algorithms.
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
void setDescription(const std::string &description)
Set algorithm description (in constructor)
const std::vector< Calibration::ExpRun > & getRunList() const
Get the list of runs for which calibration is called.
EResult
The result of calibration.
@ c_OK
Finished successfully =0 in Python.
@ c_NotEnoughData
Needs more data =2 in Python.
std::tuple< TH2F, TH2F, TH2F > DstarMassFit(std::shared_ptr< TTree > preselTree)
produce histograms for K/pi(/mu)
int m_numPBins
the number of momentum bins for the payloads
std::vector< double > CreatePBinningScheme()
build the binning scheme
bool m_isMakePlots
produce plots for monitoring
int m_MinEvtsPerTree
number of events in TTree below which we don't try to fit
int m_numDEdxBins
the number of dEdx bins for the payloads
TH2F GammaHistogram(std::shared_ptr< TTree > preselTree)
produce histograms for e
virtual EResult calibrate() override
run algorithm on data
TH2F LambdaMassFit(std::shared_ptr< TTree > preselTree)
produce histograms for protons
double m_dedxCutoff
the upper edge of the dEdx binning for the payloads
Specialized class for holding the SVD dE/dx PDFs.
Definition: SVDdEdxPDFs.h:26
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
std::vector< Atom > slice(std::vector< Atom > vec, int s, int e)
Slice the vector to contain only elements with indexes s .. e (included)
Definition: Splitter.h:85
Abstract base class for different kinds of events.