Belle II Software development
DQMHistAnalysisEcmsMonObj.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// Own header.
10#include <dqm/analysis/modules/DQMHistAnalysisEcmsMonObj.h>
11#include <framework/particledb/EvtGenDatabasePDG.h>
12#include <string>
13#include <unordered_map>
14
15#include <TROOT.h>
16#include <TRandom.h>
17#include <TH1D.h>
18#include <TLine.h>
19
20#include <RooFitResult.h>
21#include <RooRealVar.h>
22#include <RooDataSet.h>
23#include <RooGaussian.h>
24#include <RooPlot.h>
25#include <RooHist.h>
26#include <RooAddPdf.h>
27#include <RooPolyVar.h>
28#include <RooArgusBG.h>
29#include <RooSimultaneous.h>
30#include <RooCategory.h>
31#include <RooArgSet.h>
32#include <RooConstVar.h>
33
34
35using namespace std;
36using namespace Belle2;
37using namespace RooFit;
38
39//-----------------------------------------------------------------
40// Register module
41//-----------------------------------------------------------------
42
43REG_MODULE(DQMHistAnalysisEcmsMonObj);
44
47{
48 setDescription("Module to monitor Ecms information.");
49 B2DEBUG(20, "DQMHistAnalysisEcmsMonObj: Constructor done.");
50}
51
55
57{
59 B2DEBUG(20, "DQMHistAnalysisEcmsMonObj: initialized.");
60}
61
63{
64 B2DEBUG(20, "DQMHistAnalysisEcmsMonObj: beginRun called.");
65}
66
68{
69 B2DEBUG(20, "DQMHistAnalysisEcmsMonObj: event called.");
70}
71
72
73// Plot the resulting fit into the file
74TCanvas* DQMHistAnalysisEcmsMonObjModule::plotArgusFit(RooDataHist* dataE, RooAddPdf& sumB, RooArgusBG& argus, RooGaussian& gauss,
75 RooRealVar& eNow,
76 TString nTag)
77{
78 // switch to the batch mode and store the current setup
79 bool isBatch = gROOT->IsBatch();
80 gROOT->SetBatch(kTRUE);
81
82
83 // --- Plot toy data and composite PDF overlaid ---
84 RooPlot* frame = eNow.frame(40) ;
85 dataE->plotOn(frame) ;
86
87 TString name = (nTag == "B0") ? "B^{0}" : "B^{#pm}";
88 sumB.paramOn(frame, Label(name), Format("TE"), Layout(0.5, 0.8, 0.85));//, Parameters(pars));
89
90
91 sumB.plotOn(frame, Components(argus), LineStyle(kDashed)) ;
92 sumB.plotOn(frame, Components(gauss), LineStyle(kDashed), LineColor(kRed)) ;
93
94 sumB.plotOn(frame);
95
96 frame->GetXaxis()->SetTitleSize(0.0001);
97 frame->GetXaxis()->SetLabelSize(0.0001);
98 frame->SetTitleSize(0.0001);
99 frame->GetYaxis()->SetTitleSize(0.055);
100 frame->GetYaxis()->SetTitleOffset(1.1);
101 frame->GetYaxis()->SetLabelSize(0.05) ;
102
103 frame->SetTitle("");
104
105
106 RooHist* hpull = frame->pullHist() ;
107 hpull->Print();
108 RooPlot* frameRat = eNow.frame(Title("."));
109 frameRat->GetYaxis()->SetTitle("Pull");
110 frameRat->GetYaxis()->CenterTitle();
111 frameRat->GetYaxis()->SetTitleSize(0.13);
112
113 frameRat->GetYaxis()->SetNdivisions(504);
114 frameRat->GetYaxis()->SetLabelSize(0.12);
115 frameRat->GetXaxis()->SetTitleSize(0.12);
116 frameRat->GetXaxis()->SetTitleOffset(1.1);
117 frameRat->GetXaxis()->SetLabelSize(0.12);
118
119 frameRat->GetYaxis()->SetTitleOffset(0.4);
120 frameRat->addPlotable(hpull, "x0 P E1");
121 frameRat->SetMaximum(5.);
122 frameRat->SetMinimum(-5.);
123
124 TCanvas* c1 = new TCanvas(nTag + "_can");
125 TPad* pad = new TPad(nTag + "_pad", "pad", 0, 0.3, 1, 1.0);
126 pad->SetBottomMargin(0.03);
127 pad->SetRightMargin(0.05);
128 pad->SetLeftMargin(0.15);
129 pad->Draw();
130 pad->cd();
131 frame->Draw() ;
132
133
134 TLine* ll = new TLine;
135 double mRev = 10579.4e-3 / 2; // Optimal collision energy
136 ll->SetLineColor(kGreen);
137 ll->DrawLine(mRev, 0, mRev, frame->GetMaximum());
138
139
140 c1->cd();
141 TPad* padRat = new TPad(nTag + "_padRat", "padRat", 0, 0.00, 1, 0.3);
142 padRat->SetTopMargin(0.04);
143 padRat->SetBottomMargin(0.35);
144 padRat->SetRightMargin(0.05);
145 padRat->SetLeftMargin(0.15);
146 padRat->Draw();
147 padRat->cd();
148 frameRat->Draw() ;
149 c1->Update();
150
151 TLine* l = new TLine(5279.34e-3, 0.0, 5.37, 0.0);
152 l->SetLineColor(kBlue);
153 l->SetLineWidth(3);
154 l->Draw();
155
156 gROOT->SetBatch(isBatch);
157
158 return c1;
159}
160
161
162// Fit the EcmsBB histogram with Gaus+Argus function
163unordered_map<string, double> DQMHistAnalysisEcmsMonObjModule::fitEcmsBB(TH1D* hB0, TH1D* hBp)
164{
165 const double cMBp = EvtGenDatabasePDG::Instance()->GetParticle("B+")->Mass();
166 const double cMB0 = EvtGenDatabasePDG::Instance()->GetParticle("B0")->Mass();
167
168 RooRealVar eNow("eNow", "E^{*}_{B} [GeV]", cMBp, 5.37);
169
170 RooDataHist* dataE0 = new RooDataHist("dataE0hist", "Dataset from histogram", RooArgSet(eNow), hB0);
171 RooDataHist* dataEp = new RooDataHist("dataEphist", "Dataset from histogram", RooArgSet(eNow), hBp);
172
173
174 RooCategory Bcharge("sample", "sample") ;
175 Bcharge.defineType("B0") ;
176 Bcharge.defineType("Bp") ;
177
178 RooDataHist combData("combData", "combined data", eNow, Index(Bcharge), Import("B0", *dataE0), Import("Bp", *dataEp));
179
180
181 // --- Build Gaussian signal PDF ---
182 RooRealVar sigmean("mean", "#mu", 5.29, 5.27, 5.30) ;
183 RooRealVar sigwidth("sigma", "#sigma", 0.00237, 0.0001, 0.030) ;
184
185 RooGaussian gauss("gauss", "gaussian PDF", eNow, sigmean, sigwidth) ;
186
187 // --- Build Argus background PDF ---
188 RooRealVar argpar("Argus_param", "c_{Argus}", -150.7, -300., +50.0) ;
189 RooRealVar endpointBp("EndPointBp", "endPoint parameter", cMBp, 5.27, 5.291) ; //B+ value
190 RooRealVar endpointB0("EndPointB0", "endPoint parameter", cMB0, 5.27, 5.291) ; //B0 value
191 endpointB0.setConstant(kTRUE);
192 endpointBp.setConstant(kTRUE);
193
194
195 //B0 pars
196 RooPolyVar shape2B0("EndPoint2B0", "shape parameter", endpointB0, RooArgSet(RooConst(0.), RooConst(2.)));
197 RooPolyVar eNowDifB0("eNowDifB0", "eNowDifB0", eNow, RooArgSet(shape2B0, RooConst(-1.)));
198 RooArgusBG argusB0("argusB0", "Argus PDF", eNowDifB0, endpointB0, argpar) ;
199
200 //Bp pars
201 RooPolyVar shape2Bp("EndPoint2Bp", "shape parameter", endpointBp, RooArgSet(RooConst(0.), RooConst(2.)));
202 RooPolyVar eNowDifBp("eNowDifBp", "eNowDifBp", eNow, RooArgSet(shape2Bp, RooConst(-1.)));
203 RooArgusBG argusBp("argusBp", "Argus PDF", eNowDifBp, endpointBp, argpar) ;
204
205
206 // --- Construct signal+background PDF ---
207 RooRealVar nsigB0("nsigB0", "N_{sig}^{B^{0}}", 100, 0., 100000) ;
208 RooRealVar nbkgB0("nbkgB0", "N_{bg}", 100, 0., 500000) ;
209
210 RooRealVar nsigBp("nsigBp", "N_{sig}^{B^{#pm}}", 100, 0., 100000) ;
211 RooRealVar nbkgBp("nbkgBp", "N_{bg}", 100, 0., 500000) ;
212
213
214
215 RooAddPdf sumB0("sumB0", "g0+a0", RooArgList(gauss, argusB0), RooArgList(nsigB0, nbkgB0)) ;
216 RooAddPdf sumBp("sumBp", "gP+aP", RooArgList(gauss, argusBp), RooArgList(nsigBp, nbkgBp)) ;
217
218
219 // Construct a simultaneous pdf using category sample as index
220 RooSimultaneous simPdf("simPdf", "simultaneous pdf", Bcharge) ;
221
222
223 // Associate model with the physics state and model_ctl with the control state
224 simPdf.addPdf(sumB0, "B0") ;
225 simPdf.addPdf(sumBp, "Bp") ;
226
227
228 simPdf.fitTo(combData);
229
230 // Plot the results into the Canvas
231 m_canvas = new TCanvas("ecms", "ecms", 1500, 800);
232 TCanvas* c0 = plotArgusFit(dataE0, sumB0, argusB0, gauss, eNow, "B0");
233 TCanvas* cp = plotArgusFit(dataEp, sumBp, argusBp, gauss, eNow, "Bp");
234
235 m_canvas->cd();
236 m_canvas->Divide(2, 1);
237
238 m_canvas->cd(1);
239 for (auto obj : *cp->GetListOfPrimitives()) obj->DrawClone();
240 m_canvas->cd(2);
241 for (auto obj : *c0->GetListOfPrimitives()) obj->DrawClone();
242
243 // Go back to the main cBase canvas
244 m_canvas->cd();
245 m_monObj->addCanvas(m_canvas);
246
247
248 // Delete RooFit objects
249 delete dataE0;
250 delete dataEp;
251
252
253 return {{"EcmsBBcnt", 2 * sigmean.getValV()}, {"EcmsBBcntUnc", 2 * sigmean.getError()},
254 {"EcmsBBspread", 2 * sigwidth.getValV()}, {"EcmsBBspreadUnc", 2 * sigwidth.getError()}};
255}
256
257
259{
260 B2DEBUG(20, "DQMHistAnalysisEcmsMonObj: endRun called.");
261
262 auto* hB0 = (TH1D*)findHist("PhysicsObjectsMiraBelleEcmsBB/hB0");
263 auto* hBp = (TH1D*)findHist("PhysicsObjectsMiraBelleEcmsBB/hBp");
264
265 auto res = fitEcmsBB(hB0, hBp);
266
267 m_monObj->setVariable("EcmsBBcnt", res.at("EcmsBBcnt"), res.at("EcmsBBcntUnc"));
268 m_monObj->setVariable("EcmsBBspread", res.at("EcmsBBspread"), res.at("EcmsBBspreadUnc"));
269}
270
272{
273
274 B2DEBUG(20, "terminate called");
275}
276
void initialize() override final
Initialize the Module.
TCanvas * plotArgusFit(RooDataHist *dataE0, RooAddPdf &sumB0, RooArgusBG &argus, RooGaussian &gauss, RooRealVar &eNow, TString nTag="")
Plot the fit and return TCanvas with the plot.
std::unordered_map< std::string, double > fitEcmsBB(TH1D *hB0, TH1D *hBp)
Fit the histograms and return the fitted parameters.
TCanvas * m_canvas
Canvas to keep plots of the fit.
MonitoringObject * m_monObj
monitoring object
void terminate() override final
Termination action.
void event() override final
This method is called for each event.
void endRun() override final
End-of-run action.
void beginRun() override final
Called when entering a new run.
static MonitoringObject * getMonitoringObject(const std::string &name)
Get MonitoringObject with given name (new object is created if non-existing)
static TH1 * findHist(const std::string &histname, bool onlyIfUpdated=false)
Get histogram from list (no other search).
DQMHistAnalysisModule()
Constructor / Destructor.
static EvtGenDatabasePDG * Instance()
Instance method that loads the EvtGen table.
void setDescription(const std::string &description)
Sets the description of the module.
Definition Module.cc:214
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition Module.h:649
Abstract base class for different kinds of events.
STL namespace.