Belle II Software development
DQMHistAnalysisPhysics.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// File : DQMHistAnalysisPhysics.cc
10// Description : DQM module, for the physics histograms at hlt level
11//-
12
13
14#include <dqm/analysis/modules/DQMHistAnalysisPhysics.h>
15#include <TROOT.h>
16#include <RooRealVar.h>
17#include <RooDataHist.h>
18#include "RooCBShape.h"
19#include "RooCrystalBall.h"
20#include "RooChebychev.h"
21#include "RooAddPdf.h"
22#include "RooPlot.h"
23#include "RooFitResult.h"
24
25using namespace std;
26using namespace Belle2;
27
28//-----------------------------------------------------------------
29// Register the Module
30//-----------------------------------------------------------------
31REG_MODULE(DQMHistAnalysisPhysics);
32
33//-----------------------------------------------------------------
34// Implementation
35//-----------------------------------------------------------------
36
38{
39 // Description
40 // Parameter definition
41 setDescription("DQM Analysis for Physics histograms");
43 addParam("PVPrefix", m_pvPrefix, "PV Prefix", std::string("Physics:"));
44 addParam("minEntriesUPSmumu", m_minEntriesUPSmumu, "minimum number of new entries for last time slot for Ups(mumu)", 1000);
45 addParam("minEntriesUPSee", m_minEntriesUPSee, "minimum number of new entries for last time slot for Ups(ee)", 1000);
46
47 B2DEBUG(20, "DQMHistAnalysisPhysics: Constructor done.");
48}
49
53
55{
56 gROOT->cd();
57 m_monObj = getMonitoringObject("physics_hlt");
58 //new text for plots
59 m_cmUPSmumu_text = new TPaveText(0.7, 0.8, 0.9, 0.9, "NDC");
60 m_cmUPSmumu_text->SetFillColor(0);
61 m_cmUPSmumu_text->SetFillStyle(0);
62 m_cmUPSmumu_text->SetTextAlign(12);
63 m_cmUPSmumu_text->SetBorderSize(0);
64 m_cmUPSmumu_text->SetTextSize(0.026);
65 m_cmUPSee_text = new TPaveText(0.7, 0.8, 0.9, 0.9, "NDC");
66 m_cmUPSee_text->SetFillColor(0);
67 m_cmUPSee_text->SetFillStyle(0);
68 m_cmUPSee_text->SetTextAlign(12);
69 m_cmUPSee_text->SetBorderSize(0);
70 m_cmUPSee_text->SetTextSize(0.026);
71 m_ratio_text = new TPaveText(0.55, 0.6, 0.85, 0.9, "NDC");
72 m_ratio_text->SetFillColor(0);
73 m_ratio_text->SetFillStyle(0);
74 m_ratio_text->SetTextAlign(12);
75 m_ratio_text->SetBorderSize(0);
76 m_ratio_text->SetTextSize(0.026);
77
78 addDeltaPar("PhysicsObjects", "mUPSe", HistDelta::c_Entries, m_minEntriesUPSee, 1);
79 registerEpicsPV(m_pvPrefix + "mUPSee_mean", "mUPSee_mean");
80 registerEpicsPV(m_pvPrefix + "mUPSee_width", "mUPSee_width");
81 m_cmUPSee = new TCanvas("PhysicsObjects/fit_mUPSee");
82
83 addDeltaPar("PhysicsObjects", "mUPS", HistDelta::c_Entries, m_minEntriesUPSmumu, 1);
84 registerEpicsPV(m_pvPrefix + "mUPSmumu_mean", "mUPSmumu_mean");
85 registerEpicsPV(m_pvPrefix + "mUPSmumu_width", "mUPSmumu_width");
86 m_cmUPSmumu = new TCanvas("PhysicsObjects/fit_mUPSmumu");
87
88 //new ratio hadronb2_tight/bhabha_all
89 addDeltaPar("PhysicsObjects", "physicsresults", HistDelta::c_Events, 3000000, 1); // update each 10000 events
90 registerEpicsPV(m_pvPrefix + "hadronb2_tight_over_bhabha_all", "hadronb2_tight_over_bhabha_all");
91
92 B2DEBUG(20, "DQMHistAnalysisPhysics: initialized.");
93}
94
95
97{
98 B2DEBUG(20, "DQMHistAnalysisPhysics: beginRun called.");
99}
100
101void DQMHistAnalysisPhysicsModule::fitUpsilonFromHisto(TH1* histo, TPaveText* text, std::string parts, std::string prefix,
102 std::string pvname)
103{
104 double xMin = histo->GetXaxis()->GetXmin();
105 double xMax = histo->GetXaxis()->GetXmax();
106 RooRealVar m("m", parts.c_str(), xMin, xMax);
107
108 RooDataHist data("data", "histogram", m, histo);
109
110 RooRealVar mean("mean", "mass", 10.579, 10.5, 10.7);
111 RooRealVar sigma("sigma", "resolution", 0.05, 0.01, 0.15);
112 RooRealVar alphaL("alphaL", "left tail alpha", 1.5, 0.1, 5);
113 RooRealVar nL("nL", "left tail n", 1.0, 0.1, 10);
114 RooRealVar alphaR("alphaR", "right tail alpha", 2.0, 0.1, 5);
115 RooRealVar nR("nR", "right tail n", 1.0, 0.1, 10);
116
117 RooCrystalBall signal("signal", "Double CB", m, mean, sigma, sigma, alphaL, nL, alphaR, nR);
118
119 RooRealVar a0("a0", "poly constant", 0.0, -1.0, 1.0);
120 RooChebychev bkg("bkg", "Background PDF", m, RooArgList(a0));
121
122 // Set the max to 2x the total entries to ensure the fit never hits a ceiling
123 double nMax = histo->GetEntries() * 2.0;
124 RooRealVar nsig("nsig", "yield of signal", histo->GetEntries() * 0.3, 0, nMax);
125 RooRealVar nbkg("nbkg", "yield of background", histo->GetEntries() * 0.7, 0, nMax);
126
127 RooAddPdf model("model", "sig+bkg", RooArgList(signal, bkg), RooArgList(nsig, nbkg));
128
129 model.fitTo(data, RooFit::Extended(kTRUE));
130
131 RooPlot* frame = m.frame();
132 data.plotOn(frame, RooFit::DrawOption("B"), RooFit::FillColor(kGray), RooFit::LineWidth(2), RooFit::MarkerSize(0),
133 RooFit::XErrorSize(0), RooFit::DataError(RooAbsData::None), RooFit::LineStyle(kSolid), RooFit::LineColor(kBlack));
134 model.plotOn(frame, RooFit::LineColor(kBlue));
135 frame->SetTitle("Upsilon(4S) Mass Fit");
136 frame->GetXaxis()->SetTitle((parts).c_str());
137 frame->Draw();
138
139 auto measuredMass = mean.getVal();
140 auto massUncertainty = mean.getError();
141 auto measuredWidth = sigma.getVal();
142
143 text->Clear();
144 double mean_mUPS = histo->GetMean();
145 text->AddText(Form("mean : %.3f", float(mean_mUPS)));
146 text->AddText(Form("fit mean : %.3f +-%.3f", float(measuredMass), float(massUncertainty)));
147 text->AddText(Form("fit width : %.3f", float(measuredWidth)));
148 text->Draw();
149 setEpicsPV(pvname + "_mean", measuredMass);
150 setEpicsPV(pvname + "_width", measuredWidth);
151 m_monObj->setVariable(prefix + "mass", measuredMass, massUncertainty);
152 m_monObj->setVariable(prefix + "width", measuredWidth);
153
154}
155
157{
158
159
160 bool m_IsPhysicsRun = (getRunType() == "physics") || (getRunType() == "debug");
161 if (m_IsPhysicsRun == true) {
162
163 m_ratio_text->Clear();
164 auto m_hphysicsresults = findHist("PhysicsObjects/physicsresults", true);// check if updated
165 if (m_hphysicsresults) {
166 double had_ntot = m_hphysicsresults->GetBinContent(2);
167 double hadb2_ntot = m_hphysicsresults->GetBinContent(3);
168 double hadb2_tight_ntot = m_hphysicsresults->GetBinContent(4);
169 double mumu_tight_ntot = m_hphysicsresults->GetBinContent(5);
170 double bhabha_all_ntot = m_hphysicsresults->GetBinContent(6);
171
172 double ratio_hadron_bhabha = 0.;
173 double ratio_hadronb2_bhabha = 0.;
174 double ratio_hadronb2_tight_bhabha = 0.;
175 double ratio_mumu_tight_bhabha = 0.;
176 double error_hadron_bhabha = -10.;
177 double error_hadronb2_bhabha = -10.;
178 double error_hadronb2_tight_bhabha = -10.;
179 double error_mumu_tight_bhabha = -10.;
180
181
182 if (bhabha_all_ntot != 0) {
183 ratio_hadron_bhabha = had_ntot / bhabha_all_ntot;
184 error_hadron_bhabha = ratio_hadron_bhabha * sqrt((1 / had_ntot) + (1 / bhabha_all_ntot));
185 ratio_hadronb2_bhabha = hadb2_ntot / bhabha_all_ntot;
186 error_hadronb2_bhabha = ratio_hadronb2_bhabha * sqrt((1 / hadb2_ntot) + (1 / bhabha_all_ntot));
187 ratio_hadronb2_tight_bhabha = hadb2_tight_ntot / bhabha_all_ntot;
188 error_hadronb2_tight_bhabha = ratio_hadronb2_tight_bhabha * sqrt((1 / hadb2_tight_ntot) + (1 / bhabha_all_ntot));
189 ratio_mumu_tight_bhabha = mumu_tight_ntot / bhabha_all_ntot;
190 error_mumu_tight_bhabha = ratio_mumu_tight_bhabha * sqrt((1 / mumu_tight_ntot) + (1 / bhabha_all_ntot));
191 }
192 m_ratio_text->AddText(Form("hadronb2_tight/bhabha: %.4f +/- %.4f", float(ratio_hadronb2_tight_bhabha),
193 float(error_hadronb2_tight_bhabha)));
194 m_ratio_text->AddText(Form("hadronb2/bhabha: %.4f +/- %.4f", float(ratio_hadronb2_bhabha), float(error_hadronb2_bhabha)));
195 m_ratio_text->AddText(Form("mumu_tight/bhabha: %.4f +/- %.4f", float(ratio_mumu_tight_bhabha), float(error_mumu_tight_bhabha)));
196 m_ratio_text->AddText(Form("hadron/bhabha: %.4f +/- %.4f", float(ratio_hadron_bhabha), float(error_hadron_bhabha)));
197
198 }
199
200 // for pv #new hadronb2_tight/#bhabha_all
201 auto hist_hadronb2_tight_over_bhabha_all = getDelta("PhysicsObjects", "physicsresults", 0, true);// only if updated
202 if (hist_hadronb2_tight_over_bhabha_all) {
203 if (hist_hadronb2_tight_over_bhabha_all->GetBinContent(6) != 0) {
204 double hadronb2_tight_over_bhabha_all = hist_hadronb2_tight_over_bhabha_all->GetBinContent(4) /
205 hist_hadronb2_tight_over_bhabha_all->GetBinContent(6);
206 B2DEBUG(1, "hadronb2_tight_over_bhabha_all:" << hadronb2_tight_over_bhabha_all);
207 setEpicsPV("hadronb2_tight_over_bhabha_all", hadronb2_tight_over_bhabha_all);
208 }
209 }
210
211 auto hmUPSmumu = getDelta("PhysicsObjects/mUPS");
212 if (m_cmUPSmumu and hmUPSmumu) {
213 m_cmUPSmumu->cd();
214 fitUpsilonFromHisto(hmUPSmumu, m_cmUPSmumu_text, "M(#mu#mu) [GeV/c^2]", "UPSmumu", m_pvPrefix + "mUPSmumu");
215 m_cmUPSmumu->Modified();
216 m_cmUPSmumu->Update();
218 }
219 auto hmUPSee = getDelta("PhysicsObjects/mUPSe");// check if updated
220 if (m_cmUPSee and hmUPSee) {
221 m_cmUPSee->cd();
222 fitUpsilonFromHisto(hmUPSee, m_cmUPSee_text, "M(ee) [GeV/c^2]", "UPSee", m_pvPrefix + "mUPSee");
223 m_cmUPSee->Modified();
224 m_cmUPSee->Update();
226 }
227 auto* m_cphysicsresults = findCanvas("PhysicsObjects/c_physicsresults");
228 if (m_cphysicsresults) {
229 m_cphysicsresults->cd();
230 m_ratio_text->Draw();
231 m_cphysicsresults->Modified();
232 m_cphysicsresults->Update();
233 }
234 }
235}
237{
238 auto m_hphysicsresults = findHist("PhysicsObjects/physicsresults");
239 if (m_hphysicsresults) {
240 double had_ntot = m_hphysicsresults->GetBinContent(2);
241 double hadb2_ntot = m_hphysicsresults->GetBinContent(3);
242 double hadb2_tight_ntot = m_hphysicsresults->GetBinContent(4);
243 double mumu_tight_ntot = m_hphysicsresults->GetBinContent(5);
244 double bhabha_all_ntot = m_hphysicsresults->GetBinContent(6);
245 double ratio_hadron_bhabha_final = 0.;
246 double ratio_hadronb2_bhabha_final = 0.;
247 double ratio_hadronb2_tight_bhabha_final = 0.;
248 double ratio_mumu_tight_bhabha_final = 0.;
249
250 if (bhabha_all_ntot != 0) {
251 ratio_hadron_bhabha_final = had_ntot / bhabha_all_ntot;
252 ratio_hadronb2_bhabha_final = hadb2_ntot / bhabha_all_ntot;
253 ratio_hadronb2_tight_bhabha_final = hadb2_tight_ntot / bhabha_all_ntot;
254 ratio_mumu_tight_bhabha_final = mumu_tight_ntot / bhabha_all_ntot;
255 }
256 m_monObj->setVariable("ratio_hadron_bhabha_hlt", ratio_hadron_bhabha_final);
257 m_monObj->setVariable("ratio_hadronb2_bhabha_hlt", ratio_hadronb2_bhabha_final);
258 m_monObj->setVariable("ratio_hadronb2_tight_bhabha_hlt", ratio_hadronb2_tight_bhabha_final);
259 m_monObj->setVariable("ratio_mumu_tight_bhabha_hlt", ratio_mumu_tight_bhabha_final);
260 m_monObj->setVariable("hadronb2_tight_hlt", hadb2_tight_ntot);
261 m_monObj->setVariable("bhabha_all_hlt", bhabha_all_ntot);
262
263 }
264}
265
267{
268
269 B2DEBUG(20, "DQMHistAnalysisPhysics: terminate called");
270}
271
TCanvas * findCanvas(TString cname)
Find canvas by name.
static MonitoringObject * getMonitoringObject(const std::string &name)
Get MonitoringObject with given name (new object is created if non-existing)
void addDeltaPar(const std::string &dirname, const std::string &histname, HistDelta::EDeltaType t, int p, unsigned int a=1)
Add Delta histogram parameters.
static TH1 * findHist(const std::string &histname, bool onlyIfUpdated=false)
Get histogram from list (no other search).
static const std::string & getRunType(void)
Get the list of the reference histograms.
TH1 * getDelta(const std::string &fullname, int n=0, bool onlyIfUpdated=true)
Get Delta histogram.
void setEpicsPV(std::string keyname, double value)
Write value to a EPICS PV.
DQMHistAnalysisModule()
Constructor / Destructor.
int registerEpicsPV(std::string pvname, std::string keyname="")
EPICS related Functions.
void UpdateCanvas(std::string name, bool updated=true)
Mark canvas as updated (or not)
void terminate(void) override final
This method is called at the end of the event processing.
TPaveText * m_ratio_text
TPaveText, ratios physics results.
int m_minEntriesUPSee
min entries for Ups(ee)
void initialize(void) override final
Initializer.
void endRun(void) override final
End-of-run action.
std::string m_pvPrefix
prefix for EPICS PVs
void fitUpsilonFromHisto(TH1 *histo, TPaveText *text, std::string parts, std::string prefix, std::string pvname)
fit histogram with UPS mass
MonitoringObject * m_monObj
Monitoring Object.
TCanvas * m_cmUPSmumu
Tcanvas for mUPSmumu.
int m_minEntriesUPSmumu
min entries for Ups(mumu)
TPaveText * m_cmUPSee_text
TPaveText, Ups ee Invariant Mass.
TPaveText * m_cmUPSmumu_text
TPaveText, Ups Invariant Mass (mumu)
void beginRun(void) override final
Called when entering a new run.
void event(void) override final
This method is called for each event.
void setDescription(const std::string &description)
Sets the description of the module.
Definition Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition Module.h:80
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition Module.h:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition Module.h:649
double sqrt(double a)
sqrt for double
Definition beamHelpers.h:28
Abstract base class for different kinds of events.
STL namespace.