Belle II Software development
PXDRawDQMModule.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/PXDRawDQMModule.h>
10#include <pxd/dataobjects/PXDRawHit.h>
11#include <pxd/dataobjects/PXDRawAdc.h>
12#include <pxd/dataobjects/PXDDAQStatus.h>
13#include <rawdata/dataobjects/RawPXD.h>
14
15#include <vxd/geometry/GeoCache.h>
16#include <TDirectory.h>
17#include <boost/format.hpp>
18#include <string>
19
20using namespace std;
21using namespace Belle2;
22using namespace Belle2::PXD;
23using namespace Belle2::VXD;
24
25using boost::format;
26
27//-----------------------------------------------------------------
28// Register the Module
29//-----------------------------------------------------------------
30REG_MODULE(PXDRawDQM);
31
32//-----------------------------------------------------------------
33// Implementation
34//-----------------------------------------------------------------
35
37{
38 //Set module properties
39 setDescription("Monitor raw PXD");
41 addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms will be placed",
42 std::string("pxdraw"));
43
44 addParam("RawPXDsName", m_storeRawPxdrarrayName, "The name of the StoreArray of RawPXDs to be processed", string(""));
45 addParam("PXDRawHitsName", m_storeRawHitsName, "The name of the StoreArray of PXDRawHits to be processed", string(""));
46 addParam("PXDRawAdcsName", m_storeRawAdcsName, "The name of the StoreArray of PXDRawAdcs to be processed", string(""));
47}
48
50{
51 // Create a separate histogram directory and cd into it.
52 TDirectory* oldDir = gDirectory;
53 if (m_histogramDirectoryName != "") {
54 oldDir->mkdir(m_histogramDirectoryName.c_str());// do not use return value with ->cd(), its ZERO if dir already exists
55 oldDir->cd(m_histogramDirectoryName.c_str());
56 }
57
58 hrawPxdPackets = new TH1F("hrawPxdPackets", "Pxd Raw Packet Nr;Nr per Event", 16, 0, 16);
59 hrawPxdPacketSize = new TH1F("hrawPxdPacketSize", "Pxd Raw Packetsize;Words per packet", 1024, 0, 1024);
60 hrawPxdHitMapAll = new TH2F("hrawPxdHitMapAll",
61 "Pxd Raw Hit Map Overview;column+(ladder-1)*300+100;row+850*((layer-1)*2+(sensor-1))", 3700 / 50, 0, 3700, 3500 / 50, 0, 3500);
62 // ADC map not supported by DHC anymore ... deactivate filling, later remove
64 NULL;// new TH2F("hrawPxdAdcMapAll", "Pxd Raw Adc Map Overview;column+(ladder-1)*300+100;row+850*((layer-1)*2+(sensor-1))", 370/*0*/, 0, 3700, 350/*0*/, 0, 3500);
65
66 hrawPxdHitsCount = new TH1F("hrawPxdCount", "Pxd Raw Count ;Nr per Event", 8192, 0, 8192);
67
69 std::vector<VxdID> sensors = vxdGeometry.getListOfSensors();
70 for (const VxdID& avxdid : sensors) {
71 const auto& info = vxdGeometry.getSensorInfo(avxdid);
72 if (info.getType() != VXD::SensorInfoBase::PXD) continue;
73 //Only interested in PXD sensors
74 TString buff = (std::string)avxdid;
75 TString buffus = buff;
76 buffus.ReplaceAll(".", "_");
77
78 m_hrawPxdHitMap[avxdid] = new TH2F(("hrawPxdHitMap_" + buffus),
79 ("Pxd Raw Hit Map, " + buff + ";column;row"), 250,
80 0, 250, 768, 0, 768);
81 m_hrawPxdChargeMap[avxdid] = new TH2F(("hrawPxdChargeMap_" + buffus),
82 ("Pxd Raw Charge Map, " + buff + ";column;row"), 250, 0, 250, 768, 0, 768);
83 m_hrawPxdHitsCharge[avxdid] = new TH1F(("hrawPxdHitsCharge_" + buffus),
84 ("Pxd Raw Hit Charge, " + buff + ";Charge"), 256, 0, 256);
85 m_hrawPxdHitTimeWindow[avxdid] = new TH1F(("hrawPxdHitTimeWindow_" + buffus),
86 ("Pxd Raw Hit Time Window (framenr*192-gate_of_hit), " + buff + ";Time [a.u.]"), 2048, -256, 2048 - 256);
87 m_hrawPxdGateTimeWindow[avxdid] = new TH1F(("hrawPxdGateTimeWindow_" + buffus),
88 ("Pxd Raw Gate Time Window (framenr*192-triggergate_of_hit), " + buff + ";Time [a.u.]"), 2048, -256, 2048 - 256);
89 }
90 // cd back to root directory
91 oldDir->cd();
92}
93
95{
96 REG_HISTOGRAM
100 m_storeDAQEvtStats.isRequired();
101}
102
104{
105 // Just to make sure, reset all the histograms.
106 if (hrawPxdPackets) hrawPxdPackets->Reset();
111
112 for (auto& h : m_hrawPxdHitMap) if (h.second) h.second->Reset();
113 for (auto& h : m_hrawPxdChargeMap) if (h.second) h.second->Reset();
114 for (auto& h : m_hrawPxdHitsCharge) if (h.second) h.second->Reset();
115 for (auto& h : m_hrawPxdHitTimeWindow) if (h.second) h.second->Reset();
116 for (auto& h : m_hrawPxdGateTimeWindow) if (h.second) h.second->Reset();
117}
118
120{
121 hrawPxdPackets->Fill(m_storeRawPxdrarray.getEntries());
122
123 for (auto& it : m_storeRawPxdrarray) {
124 if (hrawPxdPacketSize) hrawPxdPacketSize->Fill(it.size());
125 }
126
127 if (hrawPxdHitsCount) hrawPxdHitsCount->Fill(m_storeRawHits.getEntries());
128
129 for (auto& it : m_storeRawHits) {
130
131 VxdID currentVxdId = it.getSensorID();
132 auto layer = currentVxdId.getLayerNumber();
133 auto ladder = currentVxdId.getLadderNumber();
134 auto sensor = currentVxdId.getSensorNumber();
135
136 // Get startrow and DheID from DAQEvtStats
137 const PXDDAQDHEStatus* dhe = (*m_storeDAQEvtStats).findDHE(currentVxdId);
138 if (dhe == nullptr) {
139 B2ERROR("No DHE found for SensorId: " << currentVxdId);
140 continue;
141 }
142 auto startGate = dhe->getTriggerGate();
143
144 if (hrawPxdHitMapAll) hrawPxdHitMapAll->Fill(it.getColumn() + ladder * 300 - 200,
145 100 + it.getRow() + 850 * (layer + layer + sensor - 3));
146
147 if (m_hrawPxdHitMap[currentVxdId]) m_hrawPxdHitMap[currentVxdId]->Fill(it.getColumn(), it.getRow());
148 if (m_hrawPxdChargeMap[currentVxdId]) m_hrawPxdChargeMap[currentVxdId]->Fill(it.getColumn(), it.getRow(), it.getCharge());
149 if (m_hrawPxdHitsCharge[currentVxdId]) m_hrawPxdHitsCharge[currentVxdId]->Fill(it.getCharge());
150 // Is this histogram necessary? we are folding with occupancy of sensor hits here
151 // Think about 1024*framenr-hit_row?
152 if (m_hrawPxdHitTimeWindow[currentVxdId]) m_hrawPxdHitTimeWindow[currentVxdId]->Fill(it.getFrameNr() * 192 - it.getRow() / 4);
153 if (m_hrawPxdGateTimeWindow[currentVxdId]) m_hrawPxdGateTimeWindow[currentVxdId]->Fill(it.getFrameNr() * 192 - startGate);
154 }
155
156 if (hrawPxdAdcMapAll) {
157 for (auto& it : m_storeRawAdcs) {
158 int dhh_id;
159 // calculate DHH id from Vxd Id
160 unsigned int layer, ladder, sensor;//, segment;
161 VxdID currentVxdId;
162 currentVxdId = it.getSensorID();
163 layer = currentVxdId.getLayerNumber();
164 ladder = currentVxdId.getLadderNumber();
165 sensor = currentVxdId.getSensorNumber();
166 // segment = currentVxdId.getSegmentNumber();// Frame nr? ... ignore
167 dhh_id = ((layer - 1) << 5) | ((ladder) << 1) | (sensor - 1);
168 if (dhh_id <= 0 || dhh_id >= 64) {
169 B2ERROR("SensorId (DHH ID) out of range: " << dhh_id);
170 continue;
171 }
172
173 unsigned int chip_offset;
174 chip_offset = it.getChip() * 64;
175 const unsigned char* data = &it.getData()[0];
176 for (int row = 0; row < 768; row++) {
177 for (int col = 0; col < 64; col++) {
178 hrawPxdAdcMapAll->Fill(col + chip_offset + ladder * 300 - 200, 100 + row + 850 * (layer + layer + sensor - 3), *(data++));
179 }
180 }
181 }
182 }
183}
HistoModule()
Constructor.
Definition HistoModule.h:32
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 PXD DAQ DHE Status class.
unsigned short getTriggerGate(void) const
get Trigger Gate
TH2F * hrawPxdAdcMapAll
Histogram Adc 2d hitmap (full frames only) (all pxd)
void initialize() override final
Initialize.
std::string m_storeRawAdcsName
RawAdcs StoreArray name.
TH1F * hrawPxdPacketSize
Histogram raw packet size.
StoreObjPtr< PXDDAQStatus > m_storeDAQEvtStats
Input array for DAQ Status.
StoreArray< RawPXD > m_storeRawPxdrarray
Storearray for raw data packets.
std::map< VxdID, TH2F * > m_hrawPxdChargeMap
Histogram 2d chargemap.
void defineHisto() override final
Define histograms.
std::string m_storeRawHitsName
PXDRawHits StoreArray name.
std::map< VxdID, TH1F * > m_hrawPxdGateTimeWindow
Histogram raw pixel trigger gate window.
std::string m_storeRawPxdrarrayName
RawPXD StoreArray name.
void event() override final
Event.
TH1F * hrawPxdHitsCount
Histogram pixelcount/?
TH2F * hrawPxdHitMapAll
Histogram 2d hitmap (all)
std::map< VxdID, TH1F * > m_hrawPxdHitTimeWindow
Histogram raw pixel hit "time" window.
PXDRawDQMModule()
Constructor defining the parameters.
StoreArray< PXDRawAdc > m_storeRawAdcs
Storearray for ADC from full frames.
std::string m_histogramDirectoryName
Name of the histogram directory in ROOT file.
std::map< VxdID, TH2F * > m_hrawPxdHitMap
Histogram 2d hitmap.
TH1F * hrawPxdPackets
Histogram number of raw packets.
void beginRun() override final
Begin run.
std::map< VxdID, TH1F * > m_hrawPxdHitsCharge
Histogram raw pixel charge.
StoreArray< PXDRawHit > m_storeRawHits
Storearray for raw pixels.
Class to facilitate easy access to sensor information of the VXD like coordinate transformations or p...
Definition GeoCache.h:38
const std::vector< VxdID > getListOfSensors() const
Get list of all sensors.
Definition GeoCache.cc:59
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a reference to the SensorInfo of a given SensorID.
Definition GeoCache.cc:67
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition GeoCache.cc:214
Class to uniquely identify a any structure of the PXD and SVD.
Definition VxdID.h:32
baseType getSensorNumber() const
Get the sensor id.
Definition VxdID.h:99
baseType getLadderNumber() const
Get the ladder id.
Definition VxdID.h:97
baseType getLayerNumber() const
Get the layer id.
Definition VxdID.h:95
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:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition Module.h:649
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:33
Abstract base class for different kinds of events.
STL namespace.