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_legProblem->Clear();
91 m_legProblem->AddText("ERROR!");
92 m_legProblem->AddText(Form("abs(Mean) > %3.1f ns", m_timeThreshold));
93
94 m_legWarning->Clear();
95 m_legWarning->AddText("WARNING!");
96
97 m_legNormal->Clear();
98 m_legNormal->AddText("TIME SHIFT UNDER LIMIT");
99 m_legNormal->AddText(Form("abs(Mean) < %3.1f ns", m_timeThreshold));
100
101 m_legLowStat->Clear();
102 m_legLowStat->AddText("Not enough statistics");
103
104 m_legEmpty->Clear();
105 m_legEmpty->AddText("No data recieved");
106}
107
109{
110 B2DEBUG(10, "DQMHistAnalysisSVDClustersOnTrack: event called.");
111
112 //find nEvents
113 TH1* hnEvnts = findHist("SVDExpReco/SVDDQM_nEvents", true);
114 if (hnEvnts == NULL) {
115 B2INFO("no events, nothing to do here");
116 return;
117 } else {
118 B2DEBUG(10, "SVDExpReco/SVDDQM_nEvents found");
119 }
120
121 string rtype = getRunType();
122 m_runtype = !rtype.empty() ? rtype.c_str() : "physics"; // per default
123
124 if (rtype.empty())
125 B2INFO("no run type found, put defaultwise physics");
126
127 TString tmp = hnEvnts->GetTitle();
128 Int_t pos = tmp.Last('~');
129 if (pos == -1) pos = 0;
130
131 TString runID = tmp(pos, tmp.Length() - pos);
132 B2INFO("DQMHistAnalysisSVDClustersOnTrackModule::runID = " << runID);
133 Float_t nEvents = hnEvnts->GetEntries();
134
135 // cluster time for clusters of track
136 double ratio3_6 = 0.;
137 TH1* m_h = findHist("SVDClsTrk/SVDTRK_ClusterTimeV456");
138
139 int status;
140 if (m_h != NULL) {
141
144 m_hClusterOnTrackTime_L456V.SetName("ClusterOnTrackTime_L456V");
145 m_hClusterOnTrackTime_L456V.SetTitle(Form("ClusterOnTrack Time L456V %s", runID.Data()));
146
147 if (nEvents > m_statThreshold)
149 else
150 status = lowStat;
151
155
156 } else {
157 B2INFO("Histogram SVDClsTrk/c_SVDTRK_ClusterTimeV456 from SVDDQMClustersOnTrack module not found!");
161 }
162
163 if (m_printCanvas)
164 m_cClusterOnTrackTime_L456V->Print("c_SVDClusterOnTrackTime_L456V.pdf");
165
166
167 // cluster time for clusters of track for 3 samples
168 if (m_3Samples) {
169 m_h = findHist("SVDClsTrk/SVDTRK_Cluster3TimeV456");
170
171 if (m_h != NULL) {
174 m_hClusterOnTrackTimeL456V3Samples.SetTitle(Form("ClusterOnTrack Time L456V 3 samples %s", runID.Data()));
175
176 if (nEvents > m_statThreshold)
178 else
179 status = lowStat;
180
184
185 } else {
186 B2INFO("Histogram SVDClsTrk/c_SVDTRK_Cluster3TimeV456 from SVDDQMClustersOnTrack module not found!");
190 }
191
192 if (m_printCanvas)
193 m_cClusterOnTrackTimeL456V3Samples->Print("c_SVDClusterOnTrack3Time_L456V.pdf");
194
195 ratio3_6 = m_hClusterOnTrackTimeL456V3Samples.GetEntries() / m_hClusterOnTrackTime_L456V.GetEntries();
196 }
197
198 setEpicsPV("ratio3_6", ratio3_6);
199}
200
202{
203 B2DEBUG(10, "DQMHistAnalysisSVDClustersOnTrack: endRun called");
204}
205
206
208{
209 B2DEBUG(10, "DQMHistAnalysisSVDClustersOnTrack: terminate called");
210
213}
214
216{
217 int status = good;
218 histo.GetXaxis()->SetRange(110, 190); // [-40 ns,40 ns]
219 Float_t mean_PeakInCenter = histo.GetMean(); //
220 histo.GetXaxis()->SetRange(); // back to [-150 ns,150 ns]
221 Float_t difference = 0;
222
223 if (m_runtype == "physics")
224 difference = fabs(mean_PeakInCenter - m_refMeanP);
225
226 else if (m_runtype == "cosmic")
227 difference = fabs(mean_PeakInCenter - m_refMeanC);
228
229 else {// taking cosmic limits
230 B2WARNING("Run type:" << m_runtype << "taken cosmics criteria");
231 difference = fabs(mean_PeakInCenter - m_refMeanC);
232 }
233
234 if (difference > m_timeThreshold)
235 status = error;
236 else
237 status = good;
238
239 return status;
240}
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.
void setEpicsPV(std::string keyname, double value)
Write value to a EPICS PV.
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.
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
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
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.
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
Class definition for common method.
TPaveText * m_legEmpty
plot legend, empty
TPaveText * m_legLowStat
plot legend, low stats
TPaveText * m_legWarning
plot legend, warning
void setStatusOfCanvas(int status, TCanvas *canvas, bool plotLeg=true, bool online=false)
set status of Canvas
TPaveText * m_legNormal
plot legend, normal
TPaveText * m_legProblem
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.