Belle II Software development
DAQMonitor.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 : DAQMonitor.cc
10// Description : Module to monitor raw data
11//-
12
13/* Own header. */
14#include <dqm/modules/DAQMonitor.h>
15
16/* ROOT headers. */
17#include <TDirectory.h>
18
19using namespace Belle2;
20
21REG_MODULE(DAQMonitor);
22
24{
25 setDescription("This module produces general DAQ DQM histograms.");
27}
28
30{
31 TDirectory* oldDir = gDirectory;
32 oldDir->mkdir("DAQ");
33 oldDir->cd("DAQ");
34 h_nEvt = new TH1F("Nevent", "Total Number of Events", 3, 0.0, 2.0);
35 h_pxdSize = new TH1F("PXDDataSize", "PXD Data Size;Size [kB];", 200, 0.0, 1000.0);
36 h_svdSize = new TH1F("SVDDataSize", "SVD Data Size;Size [kB];", 100, 0.0, 300.0);
37 h_cdcSize = new TH1F("CDCDataSize", "CDC Data Size;Size [kB];", 100, 0.0, 100.0);
38 h_topSize = new TH1F("TOPDataSize", "TOP Data Size;Size [kB];", 100, 0.0, 100.0);
39 h_arichSize = new TH1F("ARICHDataSize", "ARICH Data Size;Size [kB];", 100, 0.0, 40.0);
40 h_eclSize = new TH1F("ECLDataSize", "ECL Data Size;Size [kB];", 100, 0.0, 100.0);
41 h_klmSize = new TH1F("KLMDataSize", "KLM Data Size;Size [kB];", 100, 0.0, 40.0);
42 h_trgSize = new TH1F("TRGDataSize", "TRG Data Size;Size [kB];", 100, 0.0, 40.0);
43 h_hltSize = new TH1F("HLTDataSize", "HLT (Total - PXD) Data Size;Size [kB];", 200, 0.0, 1000.0);
44 h_totalSize = new TH1F("TotalDataSize", "Total (HLT + PXD) Data Size;Size [kB];", 200, 0.0, 2000.0);
45 h_runNr = new TH1F("hRunnr", "Run Number", 1, 0, 1); // define with a dummy bin, set correct range on first event
46
47 oldDir->cd();
48}
49
51{
52 REG_HISTOGRAM;
53 m_eventMetaData.isRequired();
54 m_pxdRaw.isOptional();
55 m_svdRaw.isOptional();
56 m_cdcRaw.isOptional();
57 m_topRaw.isOptional();
58 m_arichRaw.isOptional();
59 m_eclRaw.isOptional();
60 m_klmRaw.isOptional();
61 m_trgRaw.isOptional();
62}
63
65{
66 h_runNr->Reset();
67 h_nEvt->Reset();
68 h_pxdSize->Reset();
69 h_svdSize->Reset();
70 h_cdcSize->Reset();
71 h_topSize->Reset();
72 h_arichSize->Reset();
73 h_eclSize->Reset();
74 h_klmSize->Reset();
75 h_trgSize->Reset();
76 h_hltSize->Reset();
77 h_totalSize->Reset();
78}
79
81{
82 // Total number of events: just fill the histogram with 1
83 h_nEvt->Fill(1.0);
84
85 auto runNr = m_eventMetaData->getRun();
86 if (h_runNr->GetXaxis()->GetNbins() == 1) {
87 // this happens on the first event only.
88 // now we define a correct range for the histogram
89 h_runNr->SetBins(10, runNr - 5, runNr + 5); // just easy readable
90 // for the case that we have some histogram bleeding from another run
91 // we must set CanExtend or we can not merge histograms
92 // because of different axis. if th1 histogram is rebinned,
93 // we may lose the lowest bits of th1 run nr, but we still
94 // would know that we have more than one run in it!
95 // (which is the main purpose of this histogram)
96 h_runNr->GetXaxis()->SetCanExtend(kTRUE);
97 }
98 h_runNr->Fill(runNr);
99
100 // Since sizeof returns the size in bytes (B),
101 // if we divide it by 1000 we obtain kilobytes (kB).
102
103 // PXD
104 int pxdSize{0};
105 for (RawPXD& pxdRaw : m_pxdRaw)
106 pxdSize += (pxdRaw.size()) * sizeof(unsigned int);
107 h_pxdSize->Fill(static_cast<float>(pxdSize) / 1000.);
108
109 // SVD
110 int svdSize{0};
111 for (RawSVD& svdRaw : m_svdRaw) // Loop over COPPERs
112 svdSize += svdRaw.GetBlockNwords(0) * sizeof(unsigned int);
113 h_svdSize->Fill(static_cast<float>(svdSize) / 1000.);
114
115 // CDC
116 int cdcSize{0};
117 for (RawCDC& cdcRaw : m_cdcRaw) // Loop over COPPERs
118 cdcSize += cdcRaw.GetBlockNwords(0) * sizeof(unsigned int);
119 h_cdcSize->Fill(static_cast<float>(cdcSize) / 1000.);
120
121 // TOP
122 int topSize{0};
123 for (RawTOP& topRaw : m_topRaw) // Loop over COPPERs
124 topSize += topRaw.GetBlockNwords(0) * sizeof(unsigned int);
125 h_topSize->Fill(static_cast<float>(topSize) / 1000.);
126
127 // ARICH
128 int arichSize{0};
129 for (RawARICH& arichRaw : m_arichRaw) // Loop over COPPERs
130 arichSize += arichRaw.GetBlockNwords(0) * sizeof(unsigned int);
131 h_arichSize->Fill(static_cast<float>(arichSize) / 1000.);
132
133 // ECL
134 int eclSize{0};
135 for (RawECL& eclRaw : m_eclRaw) // Loop over COPPERs
136 eclSize += eclRaw.GetBlockNwords(0) * sizeof(unsigned int);
137 h_eclSize->Fill(static_cast<float>(eclSize) / 1000.);
138
139 // KLM
140 int klmSize{0};
141 for (RawKLM& klmRaw : m_klmRaw) // Loop over COPPERs
142 klmSize += klmRaw.GetBlockNwords(0) * sizeof(unsigned int);
143 h_klmSize->Fill(static_cast<float>(klmSize) / 1000.);
144
145 // TRG
146 int trgSize{0};
147 for (RawTRG& trgRaw : m_trgRaw) // Loop over COPPERs
148 trgSize += trgRaw.GetBlockNwords(0) * sizeof(unsigned int);
149 h_trgSize->Fill(static_cast<float>(trgSize) / 1000.);
150
151 // HLT size and total (HLT + PXD) size
152 // this ignores aux data (e.g. softwaretrigger result) from HLT
153 int hltSize = svdSize + cdcSize + topSize + arichSize + eclSize + klmSize + trgSize;
154 h_hltSize->Fill(static_cast<float>(hltSize) / 1000.);
155 int totalSize = pxdSize + hltSize;
156 h_totalSize->Fill(static_cast<float>(totalSize) / 1000.);
157}
TH1F * h_hltSize
Histogram for HLT data size.
Definition: DAQMonitor.h:91
void initialize() override final
Initialize.
Definition: DAQMonitor.cc:50
StoreArray< RawARICH > m_arichRaw
ARICH raw data.
Definition: DAQMonitor.h:115
DAQMonitorModule()
Constructor.
Definition: DAQMonitor.cc:23
TH1F * h_totalSize
Histogram for total data size.
Definition: DAQMonitor.h:94
StoreArray< RawSVD > m_svdRaw
SVD raw data.
Definition: DAQMonitor.h:106
StoreArray< RawTOP > m_topRaw
TOP raw data.
Definition: DAQMonitor.h:112
StoreArray< RawKLM > m_klmRaw
KLM raw data.
Definition: DAQMonitor.h:121
void defineHisto() override final
Histograms definition.
Definition: DAQMonitor.cc:29
TH1F * h_nEvt
Histogram for total number of events.
Definition: DAQMonitor.h:64
StoreArray< RawPXD > m_pxdRaw
PXD raw data.
Definition: DAQMonitor.h:103
TH1F * h_cdcSize
Histogram for CDC data size.
Definition: DAQMonitor.h:73
StoreObjPtr< EventMetaData > m_eventMetaData
Input ptr for EventMetaData.
Definition: DAQMonitor.h:100
TH1F * h_arichSize
Histogram for ARICH data size.
Definition: DAQMonitor.h:79
void event() override final
Event.
Definition: DAQMonitor.cc:80
TH1F * h_svdSize
Histogram for SVD data size.
Definition: DAQMonitor.h:70
TH1F * h_pxdSize
Histogram for PXD data size.
Definition: DAQMonitor.h:67
TH1F * h_runNr
Histogram for run nr crosscheck.
Definition: DAQMonitor.h:97
void beginRun() override final
Begin run.
Definition: DAQMonitor.cc:64
TH1F * h_trgSize
Histogram for TRG data size.
Definition: DAQMonitor.h:88
TH1F * h_eclSize
Histogram for ECL data size.
Definition: DAQMonitor.h:82
StoreArray< RawTRG > m_trgRaw
TRG raw data.
Definition: DAQMonitor.h:124
TH1F * h_klmSize
Histogram for KLM data size.
Definition: DAQMonitor.h:85
StoreArray< RawECL > m_eclRaw
ECL raw data.
Definition: DAQMonitor.h:118
StoreArray< RawCDC > m_cdcRaw
CDC raw data.
Definition: DAQMonitor.h:109
TH1F * h_topSize
Histogram for TOP data size.
Definition: DAQMonitor.h:76
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 Raw ARICH class Class for RawCOPPER class data taken by ARICH Currently, this class is almost sam...
Definition: RawARICH.h:27
The Raw CDC class Class for RawCOPPER class data taken by CDC Currently, this class is almost same as...
Definition: RawCDC.h:27
The Raw ECL class Class for RawCOPPER class data taken by ECL Currently, this class is almost same as...
Definition: RawECL.h:26
The Raw KLM class Class for RawCOPPER class data taken by KLM.
Definition: RawKLM.h:27
The Raw PXD class.
Definition: RawPXD.h:27
The Raw SVD class Class for RawCOPPER class data taken by SVD Currently, this class is almost same as...
Definition: RawSVD.h:26
The Raw TOP class Class for RawCOPPER class data taken by TOP Currently, this class is almost same as...
Definition: RawTOP.h:27
The Raw TOP class Class for RawCOPPER class data taken by TOP Currently, this class is almost same as...
Definition: RawTRG.h:27
#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.