Belle II Software development
PXDROIDQMModule.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 <pxd/modules/pxdDQM/PXDROIDQMModule.h>
10#include "TDirectory.h"
11
12using namespace std;
13using namespace Belle2;
14using namespace Belle2::PXD;
15using namespace Belle2::VXD;
16
17//-----------------------------------------------------------------
18// Register the Module
19//-----------------------------------------------------------------
20REG_MODULE(PXDROIDQM);
21
22//-----------------------------------------------------------------
23// Implementation
24//-----------------------------------------------------------------
25
26PXDROIDQMModule::PXDROIDQMModule() : HistoModule(), m_vxdGeometry(VXD::GeoCache::getInstance())
27{
28 //Set module properties
29 setDescription("Monitor ROIs");
31 addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms will be placed",
32 std::string("pxdrawroi"));
33 addParam("PXDRawROIsName", m_PXDRawROIsName, "The name of the StoreArray of PXDRawROIs to be processed", std::string(""));
34 addParam("eachModule", m_eachModule, "create for each module", false);
35 addParam("offlineDQM", m_offlineDQM, "offline DQM, use 2D plots", false);
36}
37
39{
40 TDirectory* oldDir = gDirectory;
41 if (m_histogramDirectoryName != "") {
42 oldDir->mkdir(m_histogramDirectoryName.c_str());// do not rely on return value, might be ZERO
43 oldDir->cd(m_histogramDirectoryName.c_str());//changing to the right directory
44 }
45
46 hrawROIcount = new TH1F("hrawROIcount", "ROI count;Nr per Event", 250, 0, 250);
47 hrawROItype = new TH1F("hrawROItype", "ROI type;Nr per Event", 2, 0, 2);
48
49 hrawROIHLT_DHHID = new TH1F("hrawROIHLT_DHHID", "HLT ROI per Module;DHH ID;ROIs per Module", 64, 0, 64);
50 hrawROIDC_DHHID = new TH1F("hrawROIDC_DHHID", "DATCON ROI per Module;DHH ID;ROIs per Module", 64, 0, 64);
51
52 if (m_offlineDQM) {
53 hrawROIHLTmap = new TH2F("hrawROIHLTmap", "HLT ROI Middle Map;Ucell;Vcell", 250 / 5, 0, 250, 768 / 4, 0, 768);
54 hrawROIHLTsize = new TH2F("hrawROIHLTsize", "HLT ROI Size Map;;Ucell;Vcell", 50, 0, 200, 50, 0, 200);
55 hrawROIHLTminV = new TH1F("hrawROIHLTminV", "HLT ROI minV;V", 768, 0, 768);
56 hrawROIHLTmaxV = new TH1F("hrawROIHLTmaxV", "HLT ROI maxV;V", 768, 0, 768);
57 hrawROIHLTminU = new TH1F("hrawROIHLTminU", "HLT ROI minU;U", 250, 0, 250);
58 hrawROIHLTmaxU = new TH1F("hrawROIHLTmaxU", "HLT ROI maxU;U", 250, 0, 250);
59 hrawROIHLTsizeV = new TH1F("hrawROIHLTsizeV", "HLT ROI size;V", 768, 0, 768);
60 hrawROIHLTsizeU = new TH1F("hrawROIHLTsizeU", "HLT ROI size;U", 250, 0, 250);
61
62 hrawROIDCmap = new TH2F("hrawROIDCmap", "DATCON ROI Middle Map ;U;V", 250 / 5, 0, 250, 768 / 4, 0, 768);
63 hrawROIDCsize = new TH2F("hrawROIDCsize", "DATCON ROI Size Map ;U;V", 50, 0, 200, 50, 0, 200);
64 hrawROIDCminV = new TH1F("hrawROIDCminV", "DATCON ROI minV;V", 768, 0, 768);
65 hrawROIDCmaxV = new TH1F("hrawROIDCmaxV", "DATCON ROI maxV;V", 768, 0, 768);
66 hrawROIDCminU = new TH1F("hrawROIDCminU", "DATCON ROI minU;U", 250, 0, 250);
67 hrawROIDCmaxU = new TH1F("hrawROIDCmaxU", "DATCON ROI maxU;U", 250, 0, 250);
68 hrawROIDCsizeV = new TH1F("hrawROIDCsizeV", "DATCON ROI size;V", 768, 0, 768);
69 hrawROIDCsizeU = new TH1F("hrawROIDCsizeU", "DATCON ROI size;U", 250, 0, 250);
70
71 hrawROINrDCvsNrHLT = new TH2F("hrawROINrDCvsNrHLT", "Nr DATCON ROI vs Nr HLT ROI; Nr HLT ROI;Nr DATCON ROI",
72 100, 0, 100, 100, 0, 100);
73 }
75 std::vector<VxdID> sensors = m_vxdGeometry.getListOfSensors();
76 for (VxdID& avxdid : sensors) {
78 if (info.getType() != VXD::SensorInfoBase::PXD) continue;
79 // Only interested in PXD sensors
80
81 TString buff = (std::string)avxdid;
82 TString bufful = buff;
83 bufful.ReplaceAll(".", "_");
84
85 int dhh_id = ((avxdid.getLayerNumber() - 1) << 5) | ((avxdid.getLadderNumber()) << 1) | (avxdid.getSensorNumber() - 1);
86
87 hrawROIHLTmapModule[dhh_id] = new TH2F("hrawROIHLTmap_" + bufful,
88 "HLT ROI Middle Map " + buff + ";Ucell;Vcell", 250, 0, 250, 768, 0, 768);
89 hrawROIHLTsizeModule[dhh_id] = new TH2F("hrawROIHLTsize" + bufful,
90 "HLT ROI Size Map " + buff + ";Ucell;Vcell", 50, 0, 200, 50, 0, 200);
91 hrawROIDCmapModule[dhh_id] = new TH2F("hrawROIDCmap_" + bufful,
92 "DC ROI Middle Map " + buff + ";Ucell;Vcell", 250, 0, 250, 768, 0, 768);
93 hrawROIDCsizeModule[dhh_id] = new TH2F("hrawROIDCsize" + bufful,
94 "DC ROI Size Map " + buff + ";Ucell;Vcell", 50, 0, 200, 50, 0, 200);
95 }
96 }
97
98 // cd back to root directory
99 oldDir->cd();
100}
101
103{
104 REG_HISTOGRAM
105 m_storeROIs.isRequired(m_PXDRawROIsName);
106}
107
109{
110 // Just to make sure, reset all the histograms.
111 if (hrawROIcount) hrawROIcount->Reset();
112 if (hrawROItype) hrawROItype->Reset();
114 if (hrawROIDC_DHHID) hrawROIDC_DHHID->Reset();
115
116 for (auto& a : hrawROIHLTmapModule) if (a.second) a.second->Reset();
117 for (auto& a : hrawROIHLTsizeModule) if (a.second) a.second->Reset();
118 if (hrawROIHLTmap) hrawROIHLTmap->Reset();
119 if (hrawROIHLTsize) hrawROIHLTsize->Reset();
120 if (hrawROIHLTminV) hrawROIHLTminV->Reset();
121 if (hrawROIHLTmaxV) hrawROIHLTmaxV->Reset();
122 if (hrawROIHLTminU) hrawROIHLTminU->Reset();
123 if (hrawROIHLTmaxU) hrawROIHLTmaxU->Reset();
124 if (hrawROIHLTsizeV) hrawROIHLTsizeV->Reset();
125 if (hrawROIHLTsizeU) hrawROIHLTsizeU->Reset();
126
127 for (auto& a : hrawROIDCmapModule) if (a.second) a.second->Reset();
128 for (auto& a : hrawROIDCsizeModule) if (a.second) a.second->Reset();
129 if (hrawROIDCmap) hrawROIDCmap->Reset();
130 if (hrawROIDCsize) hrawROIDCsize->Reset();
131 if (hrawROIDCminV) hrawROIDCminV->Reset();
132 if (hrawROIDCmaxV) hrawROIDCmaxV->Reset();
133 if (hrawROIDCminU) hrawROIDCminU->Reset();
134 if (hrawROIDCmaxU) hrawROIDCmaxU->Reset();
135 if (hrawROIDCsizeV) hrawROIDCsizeV->Reset();
136 if (hrawROIDCsizeU) hrawROIDCsizeU->Reset();
137
139}
140
142{
143 int nr_HLT = 0;
144 int nr_DC = 0;
145 if (hrawROIcount) hrawROIcount->Fill(-1); // misuse underflow for event count
146 for (auto& it : m_storeROIs) {
147 int nr;
148 nr = it.getNrROIs();
149 if (hrawROIcount) hrawROIcount->Fill(nr);
150 if (hrawROIDC_DHHID) hrawROIDC_DHHID->Fill(-1); // misuse underflow for roi raw paket"event" count
151 if (hrawROIHLT_DHHID) hrawROIHLT_DHHID->Fill(-1); // misuse underflow for roi raw paket"event" count
152 for (auto j = 0; j < nr; j++) {
153 if (hrawROItype) hrawROItype->Fill(it.getType(j));
154 int Vmin, Vmax, Umin, Umax, Vmean, Umean, Vsize, Usize;
155 Vmin = it.getMinVid(j);
156 Vmax = it.getMaxVid(j);
157 Vsize = Vmax - Vmin;
158 Vmean = (Vmin + Vmax) / 2;
159 Umin = it.getMinUid(j);
160 Umax = it.getMaxUid(j);
161 Usize = Umax - Umin;
162 Umean = (Umin + Umax) / 2;
163 auto id = it.getDHHID(j);
164 if (it.getType(j)) {
165 nr_DC++;
166 if (hrawROIDC_DHHID) hrawROIDC_DHHID->Fill(id);
167 if (hrawROIDCmapModule[id]) hrawROIDCmapModule[id]->Fill(Umean, Vmean);
168 if (hrawROIDCsizeModule[id]) hrawROIDCsizeModule[id]->Fill(Usize, Vsize);
169 if (hrawROIDCmap) hrawROIDCmap->Fill(Umean, Vmean);
170 if (hrawROIDCsize) hrawROIDCsize->Fill(Usize, Vsize);
171 if (hrawROIDCminV) hrawROIDCminV->Fill(Vmin);
172 if (hrawROIDCmaxV) hrawROIDCmaxV->Fill(Vmax);
173 if (hrawROIDCminU) hrawROIDCminU->Fill(Umin);
174 if (hrawROIDCmaxU) hrawROIDCmaxU->Fill(Umax);
175 if (hrawROIDCsizeV) hrawROIDCsizeV->Fill(Vsize);
176 if (hrawROIDCsizeU) hrawROIDCsizeU->Fill(Usize);
177 } else {
178 nr_HLT++;
179 if (hrawROIHLT_DHHID) hrawROIHLT_DHHID->Fill(id);
180 if (hrawROIHLTmapModule[id]) hrawROIHLTmapModule[id]->Fill(Umean, Vmean);
181 if (hrawROIHLTsizeModule[id]) hrawROIHLTsizeModule[id]->Fill(Usize, Vsize);
182 if (hrawROIHLTmap) hrawROIHLTmap->Fill(Umean, Vmean);
183 if (hrawROIHLTsize) hrawROIHLTsize->Fill(Usize, Vsize);
184 if (hrawROIHLTminV) hrawROIHLTminV->Fill(Vmin);
185 if (hrawROIHLTmaxV) hrawROIHLTmaxV->Fill(Vmax);
186 if (hrawROIHLTminU) hrawROIHLTminU->Fill(Umin);
187 if (hrawROIHLTmaxU) hrawROIHLTmaxU->Fill(Umax);
188 if (hrawROIHLTsizeV) hrawROIHLTsizeV->Fill(Vsize);
189 if (hrawROIHLTsizeU) hrawROIHLTsizeU->Fill(Usize);
190 }
191 }
192 }
193 if (hrawROINrDCvsNrHLT) hrawROINrDCvsNrHLT->Fill(nr_HLT, nr_DC);
194}
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
void initialize() override final
Initialize.
StoreArray< PXDRawROIs > m_storeROIs
Storearray for ROIs
void defineHisto() override final
Define histograms.
TH2F * hrawROINrDCvsNrHLT
Histogram
bool m_eachModule
create a histo per module
bool m_offlineDQM
create 2d plots for offline checks
void event() override final
Event.
std::string m_histogramDirectoryName
Name of the histogram directory in ROOT file.
std::map< int, TH2F * > hrawROIHLTmapModule
Histogram
void beginRun() override final
Begin run.
VXD::GeoCache & m_vxdGeometry
the VXD geometry
TH1F * hrawROIHLT_DHHID
Histogram
std::map< int, TH2F * > hrawROIHLTsizeModule
Histogram
TH1F * hrawROIcount
Histogram 2d hitmap.
PXDROIDQMModule()
Constructor defining the parameters.
std::map< int, TH2F * > hrawROIDCmapModule
Histogram
std::map< int, TH2F * > hrawROIDCsizeModule
Histogram
std::string m_PXDRawROIsName
RawROI StoreArray name.
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:39
const std::vector< VxdID > getListOfSensors() const
Get list of all sensors.
Definition: GeoCache.cc:59
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:67
Base class to provide Sensor Information for PXD and SVD.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
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
Namespace to encapsulate code needed for simulation and reconstrucion of the PXD.
Namespace to provide code needed by both Vertex Detectors, PXD and SVD, and also testbeam telescopes.
Definition: GeoCache.h:34
Abstract base class for different kinds of events.
STL namespace.