Belle II Software development
TrackingAbortDQMModule.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#include "tracking/modules/trackingDQM/TrackingAbortDQMModule.h"
10
11#include <framework/dataobjects/EventMetaData.h>
12#include <svd/dataobjects/SVDShaperDigit.h>
13#include <svd/dataobjects/SVDCluster.h>
14#include <cdc/dataobjects/CDCHit.h>
15#include <mdst/dataobjects/EventLevelTrackingInfo.h>
16#include <mdst/dataobjects/TRGSummary.h>
17
18#include <TDirectory.h>
19#include <TLine.h>
20#include <TStyle.h>
21
22#include <algorithm>
23
24
25using namespace Belle2;
26
27//-----------------------------------------------------------------
28// Register the Module
29//-----------------------------------------------------------------
30REG_MODULE(TrackingAbortDQM);
31
32
33//-----------------------------------------------------------------
34// Implementation
35//-----------------------------------------------------------------
36
38{
39 setDescription("DQM Module to monitor Tracking Aborts and detector-related quantities.");
40
41 addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms will be placed.",
42 std::string("TrackingAbort"));
43
45}
46
47
48TrackingAbortDQMModule::~TrackingAbortDQMModule()
49{
50}
51
52//------------------------------------------------------------------
53// Function to define histograms
54//-----------------------------------------------------------------
55
57{
58
59 // Create a separate histogram directories and cd into it.
60 TDirectory* oldDir = gDirectory;
61 if (m_histogramDirectoryName != "") {
62 oldDir->mkdir(m_histogramDirectoryName.c_str());
63 oldDir->cd(m_histogramDirectoryName.c_str());
64 }
65
66 //histogram index:
67 // 0 if the event is triggered OUTSIDE the active_veto window
68 std::string tag[2] = {"OUT", "IN"};
69 std::string title[2] = {"[Outside Active Veto Window]", "[Inside Active Veto Window]"};
70
71
72 //number of events with and without at least one abort
73 //outside active_veto window:
74 std::string histoName = "EventsWithAborts";
75 std::string histoTitle = "Events With at Least one Abort";
76 m_nEventsWithAbort[0] = new TH1F(TString::Format("%s_%s", histoName.c_str(), tag[0].c_str()),
77 TString::Format("%s %s", histoTitle.c_str(), title[0].c_str()),
78 2, -0.5, 1.5);
79 m_nEventsWithAbort[0]->GetYaxis()->SetTitle("Number of Events");
80 m_nEventsWithAbort[0]->GetXaxis()->SetBinLabel(1, "No Abort");
81 m_nEventsWithAbort[0]->GetXaxis()->SetBinLabel(2, "At Least One Abort");
82 m_nEventsWithAbort[0]->SetMinimum(0.1);
83
84 //inside active_veto window:
85 m_nEventsWithAbort[1] = new TH1F(*m_nEventsWithAbort[0]);
86 m_nEventsWithAbort[1]->SetName(TString::Format("%s_%s", histoName.c_str(), tag[1].c_str()));
87 m_nEventsWithAbort[1]->SetTitle(TString::Format("%s %s", histoTitle.c_str(), title[1].c_str()));
88
89 //abort flag reason
90 //outside active_veto window:
91 histoName = "TrkAbortReason";
92 histoTitle = "Tracking Abort Reason";
93 m_trackingErrorFlagsReasons[0] = new TH1F(TString::Format("%s_%s", histoName.c_str(), tag[0].c_str()),
94 TString::Format("%s %s", histoTitle.c_str(), title[0].c_str()),
95 5, -0.5, 4.5);
96 m_trackingErrorFlagsReasons[0]->GetXaxis()->SetTitle("Type of error occurred");
97 m_trackingErrorFlagsReasons[0]->GetYaxis()->SetTitle("Number of Events");
98 m_trackingErrorFlagsReasons[0]->GetXaxis()->SetBinLabel(1, "Unspecified PR");
99 m_trackingErrorFlagsReasons[0]->GetXaxis()->SetBinLabel(2, "VXDTF2");
100 m_trackingErrorFlagsReasons[0]->GetXaxis()->SetBinLabel(3, "SVDCKF");
101 m_trackingErrorFlagsReasons[0]->GetXaxis()->SetBinLabel(4, "PXDCKF");
102 m_trackingErrorFlagsReasons[0]->GetXaxis()->SetBinLabel(5, "SpacePoint");
103 //inside active_veto window:
105 m_trackingErrorFlagsReasons[1]->SetName(TString::Format("%s_%s", histoName.c_str(), tag[1].c_str()));
106 m_trackingErrorFlagsReasons[1]->SetTitle(TString::Format("%s %s", histoTitle.c_str(), title[1].c_str()));
107
108
109 //SVD L3 occupancy - see SVDDQMDose module for details
110 histoName = "SVDL3UOcc";
111 histoTitle = "SVD L3 u-side ZS5 Occupancy (%)";
112 //outside active_veto window:
113 m_svdL3uZS5Occupancy[0] = new TH1F(TString::Format("%s_%s", histoName.c_str(), tag[0].c_str()),
114 TString::Format("%s %s", histoTitle.c_str(), title[0].c_str()),
115 90, 0, 100.0 / 1536.0 * 90);
116 m_svdL3uZS5Occupancy[0]->GetXaxis()->SetTitle("occupancy [%]");
117 m_svdL3uZS5Occupancy[0]->GetYaxis()->SetTitle("Number of Events");
118 //inside active_veto window:
120 m_svdL3uZS5Occupancy[1]->SetName(TString::Format("%s_%s", histoName.c_str(), tag[1].c_str()));
121 m_svdL3uZS5Occupancy[1]->SetTitle(TString::Format("%s %s", histoTitle.c_str(), title[1].c_str()));
122
123
124 //CDC extra hits
125 histoName = "nCDCExtraHits";
126 histoTitle = "Number of CDC Extra Hits";
127 //outside active_veto window:
128 m_nCDCExtraHits[0] = new TH1F(TString::Format("%s_%s", histoName.c_str(), tag[0].c_str()),
129 TString::Format("%s %s", histoTitle.c_str(), title[0].c_str()),
130 200, 0, 5000);
131 m_nCDCExtraHits[0]->GetXaxis()->SetTitle("nCDCExtraHits");
132 m_nCDCExtraHits[0]->GetYaxis()->SetTitle("Number of Events");
133 //inside active_veto window:
134 m_nCDCExtraHits[1] = new TH1F(*m_nCDCExtraHits[0]);
135 m_nCDCExtraHits[1]->SetName(TString::Format("%s_%s", histoName.c_str(), tag[1].c_str()));
136 m_nCDCExtraHits[1]->SetTitle(TString::Format("%s %s", histoTitle.c_str(), title[1].c_str()));
137
138 //SVD L3 v-side cluster time
139 histoName = "svdL3VTime";
140 histoTitle = "Layer3 v-side Cluster Time Distribution";
141 //outside active_veto window:
142 m_svdTime[0] = new TH1F(TString::Format("%s_%s", histoName.c_str(), tag[0].c_str()),
143 TString::Format("%s %s", histoTitle.c_str(), title[0].c_str()),
144 300, -150, 150);
145 m_svdTime[0]->GetXaxis()->SetTitle("cluster time (ns)");
146 m_svdTime[0]->GetYaxis()->SetTitle("Number of Events");
147 //inside active_veto window:
148 m_svdTime[1] = new TH1F(*m_svdTime[0]);
149 m_svdTime[1]->SetName(TString::Format("%s_%s", histoName.c_str(), tag[1].c_str()));
150 m_svdTime[1]->SetTitle(TString::Format("%s %s", histoTitle.c_str(), title[1].c_str()));
151
152 //SVD, CDC Averages
153 histoName = "averages";
154 histoTitle = "Averages from SVD and CDC";
155 //outside active_veto window:
156 m_integratedAverages[0] = new TH1D(TString::Format("%s_%s", histoName.c_str(), tag[0].c_str()),
157 TString::Format("%s %s", histoTitle.c_str(), title[0].c_str()),
158 5, 0, 5);
159 m_integratedAverages[0]->GetYaxis()->SetTitle("Number of X [bin-dependent]");
160 m_integratedAverages[0]->GetXaxis()->SetBinLabel(1, "nCDCHitsInner");
161 m_integratedAverages[0]->GetXaxis()->SetBinLabel(2, "nCDCHitsOuter");
162 m_integratedAverages[0]->GetXaxis()->SetBinLabel(3, "nStripsZS5_L3V");
163 m_integratedAverages[0]->GetXaxis()->SetBinLabel(4, "nStripsZS5_L4U");
164 m_integratedAverages[0]->GetXaxis()->SetBinLabel(5, "nStripsZS5_L4V");
165 //inside active_veto window:
167 m_integratedAverages[1]->SetName(TString::Format("%s_%s", histoName.c_str(), tag[1].c_str()));
168 m_integratedAverages[1]->SetTitle(TString::Format("%s %s", histoTitle.c_str(), title[1].c_str()));
169
170 oldDir->cd();
171
172}
173
175{
176 m_eventLevelTrackingInfo.isOptional();
177 m_eventMetaData.isOptional();
178 m_trgSummary.isOptional();
179 m_strips.isOptional();
180 m_cdcHits.isOptional();
181
182 // Register histograms (calls back defineHisto)
183 REG_HISTOGRAM
184}
185
186
188{
189
190 if (m_trackingErrorFlagsReasons[0] != nullptr) m_trackingErrorFlagsReasons[0]->Reset();
191 if (m_trackingErrorFlagsReasons[1] != nullptr) m_trackingErrorFlagsReasons[1]->Reset();
192 if (m_nEventsWithAbort[0] != nullptr) m_nEventsWithAbort[0]->Reset();
193 if (m_nEventsWithAbort[1] != nullptr) m_nEventsWithAbort[1]->Reset();
194 if (m_svdL3uZS5Occupancy[0] != nullptr) m_svdL3uZS5Occupancy[0]->Reset();
195 if (m_svdL3uZS5Occupancy[1] != nullptr) m_svdL3uZS5Occupancy[1]->Reset();
196 if (m_nCDCExtraHits[0] != nullptr) m_nCDCExtraHits[0]->Reset();
197 if (m_nCDCExtraHits[1] != nullptr) m_nCDCExtraHits[1]->Reset();
198 if (m_svdTime[0] != nullptr) m_svdTime[0]->Reset();
199 if (m_svdTime[1] != nullptr) m_svdTime[1]->Reset();
200 if (m_integratedAverages[0] != nullptr) m_integratedAverages[0]->Reset();
201 if (m_integratedAverages[1] != nullptr) m_integratedAverages[1]->Reset();
202}
203
204
206{
207
208 //skip events in which we do not have EvenTMetaData or TRGSummary
209 if (!m_eventMetaData.isValid()) return;
210 if (!m_trgSummary.isValid()) return;
211
212 //skip the empty events
214 return;
215
217 return;
218
220 return;
221
223 return;
224
225 //find out if we are in the passive veto (i=0) or in the active veto window (i=1)
226 int index = 0; //events accepted in the passive veto window but not in the active
227 if (m_trgSummary->testInput("passive_veto") == 1 && m_trgSummary->testInput("cdcecl_veto") == 0) index = 1;
228
229
230 //fill the tracking abort reason histogram & nEvents with Abort
231 if (m_eventLevelTrackingInfo.isValid()) {
232 if (m_eventLevelTrackingInfo->hasAnErrorFlag()) {
233
234 m_nEventsWithAbort[index]->Fill(1);
235
236 if (m_eventLevelTrackingInfo->hasUnspecifiedTrackFindingFailure())
237 m_trackingErrorFlagsReasons[index]->Fill(0);
238 if (m_eventLevelTrackingInfo->hasVXDTF2AbortionFlag())
239 m_trackingErrorFlagsReasons[index]->Fill(1);
240 if (m_eventLevelTrackingInfo->hasSVDCKFAbortionFlag())
241 m_trackingErrorFlagsReasons[index]->Fill(2);
242 if (m_eventLevelTrackingInfo->hasPXDCKFAbortionFlag())
243 m_trackingErrorFlagsReasons[index]->Fill(3);
244 if (m_eventLevelTrackingInfo->hasSVDSpacePointCreatorAbortionFlag())
245 m_trackingErrorFlagsReasons[index]->Fill(4);
246 } else { //EventLevelTrackingInfo valid but no error
247 m_nEventsWithAbort[index]->Fill(0);
248 }
249 } else //EventLevelTrackingInfo not valid
250 m_nEventsWithAbort[index]->Fill(0);
251
252 //compute the number of ZS5 strips of L3 and L4, both sides
253 float nStripsL3UZS5 = 0;
254 float nStripsL3VZS5 = 0;
255 float nStripsL4UZS5 = 0;
256 float nStripsL4VZS5 = 0;
257 for (const SVDShaperDigit& hit : m_strips) {
258 const VxdID& sensorID = hit.getSensorID();
259 if (sensorID.getLayerNumber() > 4) continue;
260 const float noise = m_NoiseCal.getNoise(sensorID, hit.isUStrip(), hit.getCellID());
261 const float cutMinSignal = std::round(5 * noise);
262
263 if (hit.passesZS(1, cutMinSignal)) {
264 if (sensorID.getLayerNumber() == 3) {
265 if (hit.isUStrip()) nStripsL3UZS5++;
266 else nStripsL3VZS5++;
267 } else if (hit.isUStrip()) nStripsL4UZS5++;
268 else nStripsL4VZS5++;
269 }
270 }
271
272 //fill the SVD L3 v-side cluster time
273 for (const SVDCluster& hit : m_clusters) {
274 const VxdID& sensorID = hit.getSensorID();
275 if (sensorID.getLayerNumber() != 3) continue;
276 if (hit.isUCluster()) continue;
277
278 m_svdTime[index]->Fill(hit.getClsTime());
279 }
280
281 // fill the svd L3 v ZS5 occupancy, add the overflow in the last bin to make them visible in the plot
282 m_svdL3uZS5Occupancy[index]->Fill(std::min((double)nStripsL3UZS5 / m_nStripsL3U * 100, (double)5.82));
283
284 //fill the nCDCExtraHits, add the overflow in the last bin to make them visible in the plot
285 if (m_eventLevelTrackingInfo.isValid())
286 m_nCDCExtraHits[index]->Fill(std::min((int)m_eventLevelTrackingInfo->getNCDCHitsNotAssigned(), (int)4999));
287
288 //compute number of CDC hits in the inner and outer layers
289 int nCDCHitsInner = 0;
290 int nCDCHitsOuter = 0;
291 for (const CDCHit& hit : m_cdcHits) {
292 if (hit.getISuperLayer() == 0) nCDCHitsInner++;
293 else nCDCHitsOuter++;
294 }
295
296 // fill the integrated averages TH1F
297 // bin 1: nCDCHits Inner layers
298 updateBinContent(index, 1, nCDCHitsInner);
299 // bin 2: nCDCHits Outer layers
300 updateBinContent(index, 2, nCDCHitsOuter);
301 // bin 3: nStrips L3 V-side
302 updateBinContent(index, 3, nStripsL3VZS5);
303 // bin 4: nStrips L4 U-side
304 updateBinContent(index, 4, nStripsL4UZS5);
305 // bin 5: nStrips L4 V-side
306 updateBinContent(index, 5, nStripsL4VZS5);
307
308
309}
310
311void TrackingAbortDQMModule::updateBinContent(int index, int bin, float valueToBeAdded)
312{
313 float oldValue = m_integratedAverages[index]->GetBinContent(bin);
314 m_integratedAverages[index]->SetBinContent(bin, oldValue + valueToBeAdded);
315}
Class containing the result of the unpacker in raw data and the result of the digitizer in simulation...
Definition: CDCHit.h:40
@ c_B2LinkPacketCRCError
Belle2link CRC error is detected in the event.
Definition: EventMetaData.h:44
@ c_HLTCrash
The HLT reconstruction crashed in this event or the event before.
Definition: EventMetaData.h:46
@ c_ReconstructionAbort
The event was not reconstructed, e.g.
Definition: EventMetaData.h:47
@ c_B2LinkEventCRCError
HSLB_COPPER CRC error is detected in the event.
Definition: EventMetaData.h:45
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
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
The SVD Cluster class This class stores all information about reconstructed SVD clusters.
Definition: SVDCluster.h:29
float getNoise(const VxdID &sensorID, const bool &isU, const unsigned short &strip) const
This is the method for getting the noise.
The SVD ShaperDigit class.
void initialize() override final
Module function initialize.
TH1F * m_svdTime[2]
L3 V-side time for all clusters.
StoreObjPtr< EventLevelTrackingInfo > m_eventLevelTrackingInfo
tracking abort info
TH1D * m_integratedAverages[2]
integrated averages of additional SVD, CDC variables
SVDNoiseCalibrations m_NoiseCal
SVDNoise calibration db object.
TH1F * m_trackingErrorFlagsReasons[2]
stores the reason of the abort
void defineHisto() override final
Defines Histograms.
StoreObjPtr< TRGSummary > m_trgSummary
trg summary
TH1F * m_svdL3uZS5Occupancy[2]
distribution of the SVD L3 V ZS5 occupancy
StoreArray< SVDCluster > m_clusters
SVD clusters.
TH1F * m_nEventsWithAbort[2]
0: no abort; 1: at least one abort
TH1F * m_nCDCExtraHits[2]
distribution of the number of extra CDC hits
static constexpr int m_nStripsL3U
number of U-side L3 strips
StoreObjPtr< EventMetaData > m_eventMetaData
event meta data
void event() override final
Module function event.
StoreArray< SVDShaperDigit > m_strips
SVD strips.
void updateBinContent(int index, int bin, float valueToBeAdded)
function to update the bin content
std::string m_histogramDirectoryName
Name of the histogram directory in ROOT file.
StoreArray< CDCHit > m_cdcHits
CDCHits.
void beginRun() override final
Module function beginRun.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
baseType getLayerNumber() const
Get the layer id.
Definition: VxdID.h:96
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.