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
25using namespace std;
26using namespace Belle2;
27
28//-----------------------------------------------------------------
29// Register the Module
30//-----------------------------------------------------------------
31REG_MODULE(DQMHistAnalysisSVDClustersOnTrack);
32
33//-----------------------------------------------------------------
34// Implementation
35//-----------------------------------------------------------------
36
39{
40 //Parameter definition
41 B2DEBUG(10, "DQMHistAnalysisSVDClustersOnTrack: Constructor done.");
42
43 setDescription("DQM Analysis Module that produces colored canvas for a straightforward interpretation of the SVD Data Quality.");
44
45 addParam("printCanvas", m_printCanvas, "if True prints pdf of the analysis canvas", bool(false));
46 addParam("statThreshold", m_statThreshold, "Minimal number of events to compare histograms", double(10000.));
47 addParam("timeThreshold", m_timeThreshold, "Acceptable difference between mean of central peak for present and reference run",
48 double(6)); // 6 ns
49 addParam("refMCTP", m_refMeanP, "Mean of the signal time peak from Physics reference run", float(0.0)); // Approximate, from exp 20
50 addParam("refMCTC", m_refMeanC, "Mean of the signal time peak from Cosmic reference run", float(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 + "ratio3_6", "ratio3_6");
68 registerEpicsPV(m_pvPrefix + "clusterTimeOnTrackLimits", "clusTimeOnTrkLimits");
69}
70
71
73{
74 B2DEBUG(10, "DQMHistAnalysisSVDClustersOnTrack: beginRun called.");
75
77
78 if (m_3Samples)
80
81 //Retrieve limits from EPICS
82 double timeWarnUp = 0.;
83 double timeErrorLo = 0.;
84 double timeWarnLo = 0.;
85 requestLimitsFromEpicsPVs("clusTimeOnTrkLimits", timeErrorLo, timeWarnLo, timeWarnUp, m_timeThreshold);
86 B2DEBUG(10, " SVD cluster time on track threshold taken from EPICS configuration file:");
87 B2DEBUG(10, " CLUSTER TIME ON TRACK: error > " << m_timeThreshold << " ns with minimum statistics of " << m_statThreshold);
88
89 // cluster time on tracks legend
90 m_legTiProblem = new TPaveText(0.15, 0.65, 0.35, 0.80, "brNDC");
91 m_legTiProblem->AddText("ERROR!");
92 m_legTiProblem->AddText(Form("abs(Mean) > %3.1f ns", m_timeThreshold));
93 m_legTiProblem->SetFillColor(c_ColorDefault);
94 m_legTiProblem->SetLineColor(kBlack);
95
96 m_legTiNormal = new TPaveText(0.15, 0.65, 0.35, 0.80, "brNDC");
97 m_legTiNormal->AddText("TIME SHIFT UNDER LIMIT");
98 m_legTiNormal->AddText(Form("abs(Mean) < %3.1f ns", m_timeThreshold));
99 m_legTiNormal->SetFillColor(c_ColorDefault);
100 m_legTiNormal->SetLineColor(kBlack);
101
102 m_legTiEmpty = new TPaveText(0.15, 0.65, 0.35, 0.80, "brNDC");
103 m_legTiEmpty->AddText("Not enough statistics");
104 m_legTiEmpty->SetTextColor(c_ColorDefault);
105 m_legTiEmpty->SetFillColor(kBlack);
106
107 m_legTi3Problem = new TPaveText(0.15, 0.65, 0.35, 0.80, "brNDC");
108 m_legTi3Problem->AddText("ERROR!");
109 m_legTi3Problem->AddText(Form("abs(Mean) > %3.1f ns", m_timeThreshold));
110 m_legTi3Problem->SetFillColor(c_ColorDefault);
111
112 m_legTi3Normal = new TPaveText(0.15, 0.65, 0.35, 0.80, "brNDC");
113 m_legTi3Normal->AddText("TIME SHIFT UNDER LIMIT");
114 m_legTi3Normal->AddText(Form("abs(Mean) < %3.1f ns", m_timeThreshold));
115 m_legTi3Normal->SetFillColor(c_ColorDefault);
116 m_legTi3Normal->SetLineColor(kBlack);
117
118 m_legTi3Empty = new TPaveText(0.15, 0.65, 0.35, 0.80, "brNDC");
119 m_legTi3Empty->AddText("Not enough statistics");
120 m_legTi3Empty->SetTextColor(c_ColorDefault);
121 m_legTi3Empty->SetFillColor(kBlack);
122}
123
125{
126 B2DEBUG(10, "DQMHistAnalysisSVDClustersOnTrack: event called.");
127
128 //find nEvents
129 TH1* hnEvnts = findHist("SVDExpReco/SVDDQM_nEvents", true);
130 if (hnEvnts == NULL) {
131 B2INFO("no events, nothing to do here");
132 return;
133 } else {
134 B2DEBUG(10, "SVDExpReco/SVDDQM_nEvents found");
135 }
136
137 TH1* rtype = findHist("DQMInfo/rtype");
138 if (rtype)
139 B2DEBUG(10, "DQMInfo/rtype found");
140
141 m_runtype = rtype ? rtype->GetTitle() : "physics"; // per default
142
143 TString tmp = hnEvnts->GetTitle();
144 Int_t pos = tmp.Last('~');
145 if (pos == -1) pos = 0;
146
147 TString runID = tmp(pos, tmp.Length() - pos);
148 B2INFO("DQMHistAnalysisSVDClustersOnTrackModule::runID = " << runID);
149 Float_t nEvents = hnEvnts->GetEntries();
150
151 // cluster time for clusters of track
152 double ratio3_6 = 0.;
153 TH1* m_h = findHist("SVDClsTrk/SVDTRK_ClusterTimeV456");
154
155 if (m_h != NULL) {
156 bool hasError = false;
157 bool lowStat = false;
158
161 m_hClusterOnTrackTime_L456V.GetXaxis()->SetRange(110, 190); // [-40 ns,40 ns]
162 Float_t mean_PeakInCenter = m_hClusterOnTrackTime_L456V.GetMean(); //
163 m_hClusterOnTrackTime_L456V.GetXaxis()->SetRange(); // back to [-150 ns,150 ns]
164 m_hClusterOnTrackTime_L456V.SetName("ClusterOnTrackTime_L456V");
165 m_hClusterOnTrackTime_L456V.SetTitle(Form("ClusterOnTrack Time L456V %s", runID.Data()));
166
167 if (nEvents > m_statThreshold) {
168 if (m_runtype == "physics") {
169 Float_t difference_physics = fabs(mean_PeakInCenter - m_refMeanP);
170 if (difference_physics > m_timeThreshold) {
171 hasError = true;
172 }
173 } else if (m_runtype == "cosmic") {
174 Float_t difference_cosmic = fabs(mean_PeakInCenter - m_refMeanC);
175 if (difference_cosmic > m_timeThreshold) {
176 hasError = true;
177 }
178 } else { // taking cosmic limits
179 B2WARNING("Run type:" << m_runtype << "taken cosmics criteria");
180 Float_t difference_cosmic = fabs(mean_PeakInCenter - m_refMeanC);
181 if (difference_cosmic > m_timeThreshold)
182 hasError = true;
183 }
184 } else {
185 lowStat = true;
186 }
187
188 if (! hasError) {
192 m_legTiNormal->Draw();
193 } else {
197 m_legTiProblem->Draw();
198 }
199
200 if (lowStat) {
204 m_legTiEmpty->Draw();
205 }
206
207 } else {
208 B2INFO("Histogram SVDClsTrk/c_SVDTRK_ClusterTimeV456 from SVDDQMClustersOnTrack module not found!");
212 }
213
214 m_cClusterOnTrackTime_L456V->Modified();
216
217 if (m_printCanvas)
218 m_cClusterOnTrackTime_L456V->Print("c_SVDClusterOnTrackTime_L456V.pdf");
219
220
221 // cluster time for clusters of track for 3 samples
222 if (m_3Samples) {
223 m_h = findHist("SVDClsTrk/SVDTRK_Cluster3TimeV456");
224
225 if (m_h != NULL) {
226 bool hasError3 = false;
227 bool lowStat3 = false;
228
231 m_hClusterOnTrackTimeL456V3Samples.GetXaxis()->SetRange(110, 190); // [-40 ns,40 ns]
232 Float_t mean_PeakInCenter = m_hClusterOnTrackTimeL456V3Samples.GetMean(); //
233 m_hClusterOnTrackTimeL456V3Samples.GetXaxis()->SetRange(); // back to [-150 ns,150 ns]
234 m_hClusterOnTrackTimeL456V3Samples.SetTitle(Form("ClusterOnTrack Time L456V 3 samples %s", runID.Data()));
235
236 if (nEvents > m_statThreshold) {
237 if (m_runtype == "physics") {
238 Float_t difference_physics = fabs(mean_PeakInCenter - m_refMeanP);
239 if (difference_physics > m_timeThreshold) {
240 hasError3 = true;
241 }
242 } else if (m_runtype == "cosmic") {
243 Float_t difference_cosmic = fabs(mean_PeakInCenter - m_refMeanC);
244 if (difference_cosmic > m_timeThreshold) {
245 hasError3 = true;
246 }
247 } else {
248 B2WARNING("Run type:" << m_runtype);
249 }
250 } else {
251 lowStat3 = true;
252 }
253 if (! hasError3) {
257 m_legTi3Normal->Draw();
258 } else {
262 m_legTi3Problem->Draw();
263 }
264
265 if (lowStat3) {
269 m_legTi3Empty->Draw();
270 }
271
272 } else {
273 B2INFO("Histogram SVDClsTrk/c_SVDTRK_Cluster3TimeV456 from SVDDQMClustersOnTrack module not found!");
277 }
278
281
282 if (m_printCanvas)
283 m_cClusterOnTrackTimeL456V3Samples->Print("c_SVDClusterOnTrack3Time_L456V.pdf");
284
285 ratio3_6 = m_hClusterOnTrackTimeL456V3Samples.GetEntries() / m_hClusterOnTrackTime_L456V.GetEntries();
286 }
287
288 setEpicsPV("ratio3_6", ratio3_6);
289}
290
292{
293 B2DEBUG(10, "DQMHistAnalysisSVDClustersOnTrack: endRun called");
294}
295
296
298{
299 B2DEBUG(10, "DQMHistAnalysisSVDClustersOnTrack: terminate called");
300
301 delete m_legTiProblem;
302 delete m_legTiNormal;
303 delete m_legTiEmpty;
304
307}
308
The base class for the histogram analysis module.
void colorizeCanvas(TCanvas *canvas, EStatus status)
Helper function for Canvas colorization.
static TH1 * findHist(const std::string &histname, bool onlyIfUpdated=false)
Get histogram from list (no other search).
@ c_ColorDefault
default for non-coloring
void setEpicsPV(std::string keyname, double value)
Write value to a EPICS PV.
@ c_StatusDefault
default for non-coloring
@ c_StatusTooFew
Not enough entries/event to judge.
@ c_StatusError
Analysis result: Severe issue found.
@ c_StatusGood
Analysis result: Good.
int registerEpicsPV(std::string pvname, std::string keyname="")
EPICS related Functions.
bool requestLimitsFromEpicsPVs(chid id, double &lowerAlarm, double &lowerWarn, double &upperWarn, double &upperAlarm)
Get Alarm Limits from EPICS PV.
TPaveText * m_legTiEmpty
cluster time on tracks plot legend, empty
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
TPaveText * m_legTi3Normal
cluster time on tracks for 3 samples plot legend, normal
TString m_runtype
string with runtype: physics or cosmic
TPaveText * m_legTiProblem
cluster time on tracks plot legend, problem
float m_refMeanC
mean of the signal time peak from Cosmic reference run
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
TPaveText * m_legTi3Empty
cluster time on tracks for 3 samples plot legend, empty
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.
TPaveText * m_legTiNormal
cluster time on tracks plot legend, normal
float m_refMeanP
mean of the signal time peak from Physics reference run
TH1F m_hClusterOnTrackTimeL456V3Samples
time for clusters on Track for L456V histo for 3 samples
TPaveText * m_legTi3Problem
cluster time on tracks for 3 samples plot legend, problem
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:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.
STL namespace.