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
17using namespace std;
18using namespace Belle2;
19
20//-----------------------------------------------------------------
21// Register the Module
22//-----------------------------------------------------------------
23REG_MODULE(DQMHistAnalysisPhysics);
24
25//-----------------------------------------------------------------
26// Implementation
27//-----------------------------------------------------------------
28
30{
31 // Description
32 // Parameter definition
33 setDescription("DQM Analysis for Physics histograms");
35 addParam("PVPrefix", m_pvPrefix, "PV Prefix", std::string("Physics:"));
36
37 B2DEBUG(20, "DQMHistAnalysisPhysics: Constructor done.");
38}
39
41{
42}
43
45{
46 gROOT->cd();
47 m_monObj = getMonitoringObject("physics_hlt");
48 //new text for plots
49 m_cmUPS_text = new TPaveText(0.7, 0.8, 0.9, 0.9, "NDC");
50 m_cmUPS_text->SetFillColor(0);
51 m_cmUPS_text->SetFillStyle(0);
52 m_cmUPS_text->SetTextAlign(12);
53 m_cmUPS_text->SetBorderSize(0);
54 m_cmUPS_text->SetTextSize(0.026);
55 m_cmUPSe_text = new TPaveText(0.7, 0.8, 0.9, 0.9, "NDC");
56 m_cmUPSe_text->SetFillColor(0);
57 m_cmUPSe_text->SetFillStyle(0);
58 m_cmUPSe_text->SetTextAlign(12);
59 m_cmUPSe_text->SetBorderSize(0);
60 m_cmUPSe_text->SetTextSize(0.026);
61 m_ratio_text = new TPaveText(0.55, 0.6, 0.85, 0.9, "NDC");
62 m_ratio_text->SetFillColor(0);
63 m_ratio_text->SetFillStyle(0);
64 m_ratio_text->SetTextAlign(12);
65 m_ratio_text->SetBorderSize(0);
66 m_ratio_text->SetTextSize(0.026);
67
68 //new ratio hadronb2_tight/bhabha_all
69 addDeltaPar("PhysicsObjects", "physicsresults", HistDelta::c_Events, 3000000, 1); // update each 10000 events
70 registerEpicsPV(m_pvPrefix + "hadronb2_tight_over_bhabha_all", "hadronb2_tight_over_bhabha_all");
71
72 B2DEBUG(20, "DQMHistAnalysisPhysics: initialized.");
73}
74
75
77{
78 B2DEBUG(20, "DQMHistAnalysisPhysics: beginRun called.");
79}
80
82{
83
84
85 bool m_IsPhysicsRun = (getRunType() == "physics") || (getRunType() == "debug");
86 if (m_IsPhysicsRun == true) {
87
88 m_cmUPS_text->Clear();
89 auto m_hmUPS = findHist("PhysicsObjects/mUPS", true);// check if updated
90 if (m_hmUPS) {
91 double mean_mUPS = m_hmUPS->GetMean();
92 m_cmUPS_text->AddText(Form("mean : %.2f", float(mean_mUPS)));
93
94 }
95 m_cmUPSe_text->Clear();
96 auto m_hmUPSe = findHist("PhysicsObjects/mUPSe", true);// check if updated
97 if (m_hmUPSe) {
98 double mean_mUPSe = m_hmUPSe->GetMean();
99 m_cmUPSe_text->AddText(Form("mean : %.2f", float(mean_mUPSe)));
100 }
101 m_ratio_text->Clear();
102 auto m_hphysicsresults = findHist("PhysicsObjects/physicsresults", true);// check if updated
103 if (m_hphysicsresults) {
104 double had_ntot = m_hphysicsresults->GetBinContent(2);
105 double hadb2_ntot = m_hphysicsresults->GetBinContent(3);
106 double hadb2_tight_ntot = m_hphysicsresults->GetBinContent(4);
107 double mumu_tight_ntot = m_hphysicsresults->GetBinContent(5);
108 double bhabha_all_ntot = m_hphysicsresults->GetBinContent(6);
109
110 double ratio_hadron_bhabha = 0.;
111 double ratio_hadronb2_bhabha = 0.;
112 double ratio_hadronb2_tight_bhabha = 0.;
113 double ratio_mumu_tight_bhabha = 0.;
114 double error_hadron_bhabha = -10.;
115 double error_hadronb2_bhabha = -10.;
116 double error_hadronb2_tight_bhabha = -10.;
117 double error_mumu_tight_bhabha = -10.;
118
119
120 if (bhabha_all_ntot != 0) {
121 ratio_hadron_bhabha = had_ntot / bhabha_all_ntot;
122 error_hadron_bhabha = ratio_hadron_bhabha * sqrt((1 / had_ntot) + (1 / bhabha_all_ntot));
123 ratio_hadronb2_bhabha = hadb2_ntot / bhabha_all_ntot;
124 error_hadronb2_bhabha = ratio_hadronb2_bhabha * sqrt((1 / hadb2_ntot) + (1 / bhabha_all_ntot));
125 ratio_hadronb2_tight_bhabha = hadb2_tight_ntot / bhabha_all_ntot;
126 error_hadronb2_tight_bhabha = ratio_hadronb2_tight_bhabha * sqrt((1 / hadb2_tight_ntot) + (1 / bhabha_all_ntot));
127 ratio_mumu_tight_bhabha = mumu_tight_ntot / bhabha_all_ntot;
128 error_mumu_tight_bhabha = ratio_mumu_tight_bhabha * sqrt((1 / mumu_tight_ntot) + (1 / bhabha_all_ntot));
129 }
130 m_ratio_text->AddText(Form("hadronb2_tight/bhabha: %.4f +/- %.4f", float(ratio_hadronb2_tight_bhabha),
131 float(error_hadronb2_tight_bhabha)));
132 m_ratio_text->AddText(Form("hadronb2/bhabha: %.4f +/- %.4f", float(ratio_hadronb2_bhabha), float(error_hadronb2_bhabha)));
133 m_ratio_text->AddText(Form("mumu_tight/bhabha: %.4f +/- %.4f", float(ratio_mumu_tight_bhabha), float(error_mumu_tight_bhabha)));
134 m_ratio_text->AddText(Form("hadron/bhabha: %.4f +/- %.4f", float(ratio_hadron_bhabha), float(error_hadron_bhabha)));
135
136 }
137
138 // for pv #new hadronb2_tight/#bhabha_all
139 auto hist_hadronb2_tight_over_bhabha_all = getDelta("PhysicsObjects", "physicsresults", 0, true);// only if updated
140 if (hist_hadronb2_tight_over_bhabha_all) {
141 if (hist_hadronb2_tight_over_bhabha_all->GetBinContent(6) != 0) {
142 double hadronb2_tight_over_bhabha_all = hist_hadronb2_tight_over_bhabha_all->GetBinContent(4) /
143 hist_hadronb2_tight_over_bhabha_all->GetBinContent(6);
144 B2DEBUG(1, "hadronb2_tight_over_bhabha_all:" << hadronb2_tight_over_bhabha_all);
145 setEpicsPV("hadronb2_tight_over_bhabha_all", hadronb2_tight_over_bhabha_all);
146 }
147 }
148
149 auto* m_cmUPS = findCanvas("PhysicsObjects/c_mUPS");
150 if (m_cmUPS) {
151 m_cmUPS->cd();
152 m_cmUPS_text->Draw();
153 m_cmUPS->Modified();
154 m_cmUPS->Update();
155 }
156 auto* m_cmUPSe = findCanvas("PhysicsObjects/c_mUPSe");
157 if (m_cmUPSe) {
158 m_cmUPSe->cd();
159 m_cmUPSe_text->Draw();
160 m_cmUPSe->Modified();
161 m_cmUPSe->Update();
162 }
163 auto* m_cphysicsresults = findCanvas("PhysicsObjects/c_physicsresults");
164 if (m_cphysicsresults) {
165 m_cphysicsresults->cd();
166 m_ratio_text->Draw();
167 m_cphysicsresults->Modified();
168 m_cphysicsresults->Update();
169 }
170 }
171}
173{
174 auto m_hphysicsresults = findHist("PhysicsObjects/physicsresults");
175 if (m_hphysicsresults) {
176 double had_ntot = m_hphysicsresults->GetBinContent(2);
177 double hadb2_ntot = m_hphysicsresults->GetBinContent(3);
178 double hadb2_tight_ntot = m_hphysicsresults->GetBinContent(4);
179 double mumu_tight_ntot = m_hphysicsresults->GetBinContent(5);
180 double bhabha_all_ntot = m_hphysicsresults->GetBinContent(6);
181 double ratio_hadron_bhabha_final = 0.;
182 double ratio_hadronb2_bhabha_final = 0.;
183 double ratio_hadronb2_tight_bhabha_final = 0.;
184 double ratio_mumu_tight_bhabha_final = 0.;
185
186 if (bhabha_all_ntot != 0) {
187 ratio_hadron_bhabha_final = had_ntot / bhabha_all_ntot;
188 ratio_hadronb2_bhabha_final = hadb2_ntot / bhabha_all_ntot;
189 ratio_hadronb2_tight_bhabha_final = hadb2_tight_ntot / bhabha_all_ntot;
190 ratio_mumu_tight_bhabha_final = mumu_tight_ntot / bhabha_all_ntot;
191 }
192 m_monObj->setVariable("ratio_hadron_bhabha_hlt", ratio_hadron_bhabha_final);
193 m_monObj->setVariable("ratio_hadronb2_bhabha_hlt", ratio_hadronb2_bhabha_final);
194 m_monObj->setVariable("ratio_hadronb2_tight_bhabha_hlt", ratio_hadronb2_tight_bhabha_final);
195 m_monObj->setVariable("ratio_mumu_tight_bhabha_hlt", ratio_mumu_tight_bhabha_final);
196
197 }
198}
199
201{
202
203 B2DEBUG(20, "DQMHistAnalysisPhysics: terminate called");
204}
205
The base class for the histogram analysis module.
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 Run Type.
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.
int registerEpicsPV(std::string pvname, std::string keyname="")
EPICS related Functions.
void terminate(void) override final
This method is called at the end of the event processing.
TPaveText * m_ratio_text
TPaveText, ratios physics results.
void initialize(void) override final
Initializer.
void endRun(void) override final
End-of-run action.
std::string m_pvPrefix
prefix for EPICS PVs
MonitoringObject * m_monObj
Monitoring Object.
TPaveText * m_cmUPS_text
TPaveText, Ups Invariant Mass (mumu)
TPaveText * m_cmUPSe_text
TPaveText, Ups ee Invariant Mass.
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 setVariable(const std::string &var, float val, float upErr=-1., float dwErr=-1)
set value to float variable (new variable is made if not yet existing)
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.