Belle II Software development
DQMHistAnalysisSVDClustersOnTrack.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 : DQMHistAnalysisSVDClustersOnTrack.cc
10// Description : module for DQM histogram analysis of SVD sensors occupancies
11//-
12
13
14#include <dqm/analysis/modules/DQMHistAnalysisSVDClustersOnTrack.h>
15#include <vxd/geometry/GeoCache.h>
16
17#include <TROOT.h>
18#include <TStyle.h>
19#include <TString.h>
20#include <TAxis.h>
21
22#include <TMath.h>
23#include <iostream>
24#include <cmath>
25
26using namespace std;
27using namespace Belle2;
28
29//-----------------------------------------------------------------
30// Register the Module
31//-----------------------------------------------------------------
32REG_MODULE(DQMHistAnalysisSVDClustersOnTrack);
33
34//-----------------------------------------------------------------
35// Implementation
36//-----------------------------------------------------------------
37
40{
41 //Parameter definition
42 B2DEBUG(10, "DQMHistAnalysisSVDClustersOnTrack: Constructor done.");
43
44 setDescription("DQM Analysis Module that produces colored canvas for a straightforward interpretation of the SVD Data Quality.");
45
46 addParam("printCanvas", m_printCanvas, "if True prints pdf of the analysis canvas", bool(false));
47 addParam("statThreshold", m_statThreshold, "Minimal number of events to compare histograms", double(10000.));
48 addParam("timeThreshold", m_timeThreshold, "Acceptable difference between mean of central peak for present and reference run",
49 double(8)); // 8 ns
50 addParam("refMode", m_refMode, "Mode reference of the signal time peak", double(0.0));
51 addParam("samples3", m_3Samples, "if True 3 samples histograms analysis is performed", bool(false));
52 addParam("PVPrefix", m_pvPrefix, "PV Prefix", std::string("SVD:"));
53}
54
56
58{
59 B2DEBUG(10, "DQMHistAnalysisSVDClustersOnTrack: initialized.");
60
61 m_cClusterOnTrackTime_L456V = new TCanvas("SVDAnalysis/c_ClusterOnTrackTime_L456V");
62
63 if (m_3Samples)
64 m_cClusterOnTrackTimeL456V3Samples = new TCanvas("SVDAnalysis/c_ClusterOnTrackTime_L456V3Samples");
65
66 //register limits for EPICS
67 registerEpicsPV(m_pvPrefix + "clusterTimeOnTrackMode", "clusterTimeOnTrackMode");
68 registerEpicsPV(m_pvPrefix + "ratio3_6", "ratio3_6");
69 registerEpicsPV(m_pvPrefix + "clusterTimeOnTrackLimits", "clusTimeOnTrkLimits");
70}
71
72
74{
75 B2DEBUG(10, "DQMHistAnalysisSVDClustersOnTrack: beginRun called.");
76
78
79 if (m_3Samples)
81
82 //Retrieve limits from EPICS
83 double timeWarnUp = 0.;
84 double timeErrorLo = 0.;
85 double timeWarnLo = 0.;
86 requestLimitsFromEpicsPVs("clusTimeOnTrkLimits", timeErrorLo, timeWarnLo, timeWarnUp, m_timeThreshold);
87
88 //Retrieve mode value for cluster time on track
89 double ref = getEpicsPV("clusterTimeOnTrackMode");
90 if (!std::isnan(ref))
91 m_refMode = ref;
92
93 B2DEBUG(10, " SVD cluster time on track threshold taken from EPICS configuration file:");
94 B2DEBUG(10, " CLUSTER TIME ON TRACK: error > " << m_timeThreshold - m_refMode << " ns with minimum statistics of " <<
96
97 // cluster time on tracks legend
98 m_legProblem->Clear();
99 m_legProblem->AddText("ERROR!");
100 m_legProblem->AddText(Form("abs(Mode - Ref) > %3.1f ns", m_timeThreshold));
101 m_legProblem->AddText("Mode: 0.0 ns");
102
103
104 m_legWarning->Clear();
105 m_legWarning->AddText("WARNING!");
106
107 m_legNormal->Clear();
108 m_legNormal->AddText("TIME SHIFT UNDER LIMIT");
109 m_legNormal->AddText(Form("abs(Mode - Ref) < %3.1f ns", m_timeThreshold));
110 m_legNormal->AddText("Mode: 0.0 ns");
111
112 m_legLowStat->Clear();
113 m_legLowStat->AddText("Not enough statistics");
114
115 m_legEmpty->Clear();
116 m_legEmpty->AddText("No data recieved");
117}
118
120{
121 B2DEBUG(10, "DQMHistAnalysisSVDClustersOnTrack: event called.");
122
123 //find nEvents
124 TH1* hnEvnts = findHist("SVDExpReco/SVDDQM_nEvents", true);
125 if (hnEvnts == NULL) {
126 B2INFO("no events, nothing to do here");
127 return;
128 } else {
129 B2DEBUG(10, "SVDExpReco/SVDDQM_nEvents found");
130 }
131
132 string rtype = getRunType();
133 m_runtype = !rtype.empty() ? rtype.c_str() : "physics"; // per default
134
135 if (rtype.empty())
136 B2INFO("no run type found, put defaultwise physics");
137
138 TString tmp = hnEvnts->GetTitle();
139 Int_t pos = tmp.Last('~');
140 if (pos == -1) pos = 0;
141
142 TString runID = tmp(pos, tmp.Length() - pos);
143 B2INFO("DQMHistAnalysisSVDClustersOnTrackModule::runID = " << runID);
144 Float_t nEvents = hnEvnts->GetEntries();
145
146 // cluster time for clusters of track
147 double ratio3_6 = 0.;
148 TH1* m_h = findHist("SVDClsTrk/SVDTRK_ClusterTimeV456");
149
150 int status;
151 if (m_h != NULL) {
152
155 TString hName = getHistoNameFromCanvas(m_cClusterOnTrackTime_L456V->GetName());
156 m_hClusterOnTrackTime_L456V.SetName(hName.Data());
157 m_hClusterOnTrackTime_L456V.SetTitle(Form("ClusterOnTrack Time L456V %s", runID.Data()));
158 m_hClusterOnTrackTime_L456V.SetStats(false);
159
160 Int_t binMax = m_hClusterOnTrackTime_L456V.GetMaximumBin();
161 double mode = m_hClusterOnTrackTime_L456V.GetXaxis()->GetBinCenter(binMax);
162
163 if (nEvents > m_statThreshold)
164 status = getCanvasStatus(mode);
165 else
166 status = lowStat;
167
171
172 } else {
173 B2INFO("Histogram SVDClsTrk/c_SVDTRK_ClusterTimeV456 from SVDDQMClustersOnTrack module not found!");
177 }
178
179 if (m_printCanvas)
180 m_cClusterOnTrackTime_L456V->Print("c_SVDClusterOnTrackTime_L456V.pdf");
181
182
183 // cluster time for clusters of track for 3 samples
184 if (m_3Samples) {
185 m_h = findHist("SVDClsTrk/SVDTRK_Cluster3TimeV456");
186
187 if (m_h != NULL) {
191 m_hClusterOnTrackTimeL456V3Samples.SetName(hName.Data());
192 m_hClusterOnTrackTimeL456V3Samples.SetTitle(Form("ClusterOnTrack Time L456V 3 samples %s", runID.Data()));
193
194 Int_t binMax = m_hClusterOnTrackTime_L456V.GetMaximumBin();
195 double mode = m_hClusterOnTrackTime_L456V.GetXaxis()->GetBinCenter(binMax);
196
197
198 if (nEvents > m_statThreshold)
199 status = getCanvasStatus(mode);
200 else
201 status = lowStat;
202
206
207 } else {
208 B2INFO("Histogram SVDClsTrk/c_SVDTRK_Cluster3TimeV456 from SVDDQMClustersOnTrack module not found!");
212 }
213
214 if (m_printCanvas)
215 m_cClusterOnTrackTimeL456V3Samples->Print("c_SVDClusterOnTrack3Time_L456V.pdf");
216
217 ratio3_6 = m_hClusterOnTrackTimeL456V3Samples.GetEntries() / m_hClusterOnTrackTime_L456V.GetEntries();
218 }
219
220 setEpicsPV("ratio3_6", ratio3_6);
221}
222
224{
225 B2DEBUG(10, "DQMHistAnalysisSVDClustersOnTrack: endRun called");
226}
227
228
230{
231 B2DEBUG(10, "DQMHistAnalysisSVDClustersOnTrack: terminate called");
232
235}
236
238{
239 int status = good;
240
241 if (fabs(mode - m_refMode) > m_timeThreshold) {
242 status = error;
243 TText* text = m_legProblem->GetLine(m_legProblem->GetSize() - 1);
244 text->SetText(text->GetX(), text->GetY(), Form("Mode - Ref: %3.1f ns", mode - m_refMode));
245 } else {
246 status = good;
247 TText* text = m_legNormal->GetLine(m_legNormal->GetSize() - 1);
248 text->SetText(text->GetX(), text->GetY(), Form("Mode - Ref: %3.1f ns", mode - m_refMode));
249 }
250
251 return status;
252}
int registerEpicsPV(const std::string &pvname, const std::string &keyname="")
EPICS related Functions.
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.
double getEpicsPV(const std::string &keyname)
Read value from a EPICS PV.
bool requestLimitsFromEpicsPVs(chid id, double &lowerAlarm, double &lowerWarn, double &upperWarn, double &upperAlarm)
Get Alarm Limits from EPICS PV.
void setEpicsPV(const std::string &keyname, double value)
Write value to a EPICS PV.
TCanvas * m_cClusterOnTrackTimeL456V3Samples
time for clusters on Track for L456V canvas for 3 sampples
double m_statThreshold
minimal number of events to compare histograms
double m_timeThreshold
difference between mean of cluster time for present and reference run
TString m_runtype
string with runtype: physics or cosmic
double m_refMode
reference mode of the signal time peak
void terminate() override final
This method is called at the end of the event processing.
void event() override final
This method is called for each event.
bool m_printCanvas
Parameters accesible from basf2 scripts.
bool m_3Samples
if true enable 3 samples histograms analysis
TH1F m_hClusterOnTrackTime_L456V
time for clusters on Track for L456V histo
void endRun() override final
This method is called if the current run ends.
TCanvas * m_cClusterOnTrackTime_L456V
time for clusters on Track for L456V canvas
void beginRun() override final
Called when entering a new run.
TH1F m_hClusterOnTrackTimeL456V3Samples
time for clusters on Track for L456V histo for 3 samples
static TString getHistoNameFromCanvas(TString cName, TString view="", TString cPrefix="c_", TString hPrefix="")
get histogram name from Canvas name
TPaveText * m_legEmpty
plot legend, empty
TPaveText * m_legLowStat
plot legend, low stats
TPaveText * m_legWarning
plot legend, warning
TPaveText * m_legNormal
plot legend, normal
TPaveText * m_legProblem
plot legend, problem
DQMHistAnalysisSVDModule(bool panelTop=false, bool online=false, bool groupIDs=false)
Constructor.
void setStatusOfCanvas(int status, TCanvas *canvas, bool plotLeg=true, int histoType=kOffline)
set status of Canvas
void setDescription(const std::string &description)
Sets the description of the module.
Definition Module.cc:214
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
Abstract base class for different kinds of events.
STL namespace.