Belle II Software  release-08-01-10
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 <vxd/geometry/GeoCache.h>
11 
12 #include <TDirectory.h>
13 #include <boost/format.hpp>
14 #include <string>
15 
16 using namespace std;
17 using namespace Belle2;
18 using namespace Belle2::PXD;
19 using namespace Belle2::VXD;
20 
21 using boost::format;
22 
23 //-----------------------------------------------------------------
24 // Register the Module
25 //-----------------------------------------------------------------
26 REG_MODULE(PXDRawDQM);
27 
28 //-----------------------------------------------------------------
29 // Implementation
30 //-----------------------------------------------------------------
31 
32 PXDRawDQMModule::PXDRawDQMModule() : HistoModule(), m_storeRawPxdrarray(), m_storeRawHits(), m_storeRawAdcs()
33 {
34  //Set module properties
35  setDescription("Monitor raw PXD");
37  addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms will be placed",
38  std::string("pxdraw"));
39 
40  addParam("RawPXDsName", m_storeRawPxdrarrayName, "The name of the StoreArray of RawPXDs to be processed", string(""));
41  addParam("PXDRawHitsName", m_storeRawHitsName, "The name of the StoreArray of PXDRawHits to be processed", string(""));
42  addParam("PXDRawAdcsName", m_storeRawAdcsName, "The name of the StoreArray of PXDRawAdcs to be processed", string(""));
43 }
44 
46 {
47  // Create a separate histogram directory and cd into it.
48  TDirectory* oldDir = gDirectory;
49  if (m_histogramDirectoryName != "") {
50  oldDir->mkdir(m_histogramDirectoryName.c_str());// do not use return value with ->cd(), its ZERO if dir already exists
51  oldDir->cd(m_histogramDirectoryName.c_str());
52  }
53 
54  hrawPxdPackets = new TH1F("hrawPxdPackets", "Pxd Raw Packet Nr;Nr per Event", 16, 0, 16);
55  hrawPxdPacketSize = new TH1F("hrawPxdPacketSize", "Pxd Raw Packetsize;Words per packet", 1024, 0, 1024);
56  hrawPxdHitMapAll = new TH2F("hrawPxdHitMapAll",
57  "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);
58  // ADC map not supported by DHC anymore ... deactive filling, later remove
60  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);
61 
62  hrawPxdHitsCount = new TH1F("hrawPxdCount", "Pxd Raw Count ;Nr per Event", 8192, 0, 8192);
63  for (auto i = 0; i < 64; i++) {
64  auto layer = (((i >> 5) & 0x1) + 1);
65  auto ladder = ((i >> 1) & 0xF);
66  auto sensor = ((i & 0x1) + 1);
67 
68  // Check if sensor exist
69  if (Belle2::VXD::GeoCache::getInstance().validSensorID(Belle2::VxdID(layer, ladder, sensor))) {
70  string s = str(format("Sensor %d:%d:%d (DHH ID %02Xh)") % layer % ladder % sensor % i);
71  string s2 = str(format("_%d.%d.%d") % layer % ladder % sensor);
72 
73  hrawPxdHitMap[i] = new TH2F(("hrawPxdHitMap" + s2).c_str(),
74  ("Pxd Raw Hit Map, " + s + ";column;row").c_str(), 250,
75  0, 250, 768, 0, 768);
76  hrawPxdChargeMap[i] = new TH2F(("hrawPxdChargeMap" + s2).c_str(),
77  ("Pxd Raw Charge Map, " + s + ";column;row").c_str(), 250, 0, 250, 768, 0, 768);
78  hrawPxdHitsCharge[i] = new TH1F(("hrawPxdHitsCharge" + s2).c_str(),
79  ("Pxd Raw Hit Charge, " + s + ";Charge").c_str(), 256, 0, 256);
80  hrawPxdHitTimeWindow[i] = new TH1F(("hrawPxdHitTimeWindow" + s2).c_str(),
81  ("Pxd Raw Hit Time Window (framenr*192-gate_of_hit), " + s + ";Time [a.u.]").c_str(), 2048, -256, 2048 - 256);
82  hrawPxdGateTimeWindow[i] = new TH1F(("hrawPxdGateTimeWindow" + s2).c_str(),
83  ("Pxd Raw Gate Time Window (framenr*192-triggergate_of_hit), " + s + ";Time [a.u.]").c_str(), 2048, -256, 2048 - 256);
84  } else {
85  hrawPxdHitMap[i] = NULL;
86  hrawPxdChargeMap[i] = NULL;
87  hrawPxdHitsCharge[i] = NULL;
88  hrawPxdHitTimeWindow[i] = NULL;
89  hrawPxdGateTimeWindow[i] = NULL;
90  }
91  }
92 
93  // cd back to root directory
94  oldDir->cd();
95 }
96 
98 {
99  REG_HISTOGRAM
103  m_storeDAQEvtStats.isRequired();
104 }
105 
107 {
108  // Just to make sure, reset all the histograms.
109  if (hrawPxdPackets) hrawPxdPackets->Reset();
110  if (hrawPxdPacketSize) hrawPxdPacketSize->Reset();
111  if (hrawPxdHitsCount) hrawPxdHitsCount->Reset();
112  if (hrawPxdHitMapAll) hrawPxdHitMapAll->Reset();
113  if (hrawPxdAdcMapAll) hrawPxdAdcMapAll->Reset();
114  for (int i = 0; i < 64; i++) {
115  if (hrawPxdHitMap[i]) hrawPxdHitMap[i]->Reset();
116  if (hrawPxdChargeMap[i]) hrawPxdChargeMap[i]->Reset();
117  if (hrawPxdHitsCharge[i]) hrawPxdHitsCharge[i]->Reset();
118  if (hrawPxdHitTimeWindow[i]) hrawPxdHitTimeWindow[i]->Reset();
119  if (hrawPxdGateTimeWindow[i]) hrawPxdGateTimeWindow[i]->Reset();
120  }
121 }
122 
124 {
125  hrawPxdPackets->Fill(m_storeRawPxdrarray.getEntries());
126 
127  for (auto& it : m_storeRawPxdrarray) {
128  if (hrawPxdPacketSize) hrawPxdPacketSize->Fill(it.size());
129  }
130 
131  if (hrawPxdHitsCount) hrawPxdHitsCount->Fill(m_storeRawHits.getEntries());
132 
133  for (auto& it : m_storeRawHits) {
134 
135  VxdID currentVxdId = it.getSensorID();
136  auto layer = currentVxdId.getLayerNumber();
137  auto ladder = currentVxdId.getLadderNumber();
138  auto sensor = currentVxdId.getSensorNumber();
139 
140  // Get startrow and DheID from DAQEvtStats
141  const PXDDAQDHEStatus* dhe = (*m_storeDAQEvtStats).findDHE(currentVxdId);
142  if (dhe == nullptr) {
143  B2ERROR("No DHE found for SensorId: " << currentVxdId);
144  continue;
145  }
146 
147  auto dhh_id = dhe->getDHEID();
148  auto startGate = dhe->getTriggerGate();
149 
150  if (dhh_id == 0 || dhh_id >= 64) {
151  B2ERROR("SensorId (DHH ID) out of range: " << dhh_id);
152  continue;
153  }
154  if (hrawPxdHitMap[dhh_id]) hrawPxdHitMap[dhh_id]->Fill(it.getColumn(), it.getRow());
155  if (hrawPxdHitMapAll) hrawPxdHitMapAll->Fill(it.getColumn() + ladder * 300 - 200,
156  100 + it.getRow() + 850 * (layer + layer + sensor - 3));
157  if (hrawPxdChargeMap[dhh_id]) hrawPxdChargeMap[dhh_id]->Fill(it.getColumn(), it.getRow(), it.getCharge());
158  if (hrawPxdHitsCharge[dhh_id]) hrawPxdHitsCharge[dhh_id]->Fill(it.getCharge());
159  // Is this histogram necessary? we are folding with occupancy of sensor hits here
160  // Think about 1024*framenr-hit_row?
161  if (hrawPxdHitTimeWindow[dhh_id]) hrawPxdHitTimeWindow[dhh_id]->Fill(it.getFrameNr() * 192 - it.getRow() / 4);
162  if (hrawPxdGateTimeWindow[dhh_id]) hrawPxdGateTimeWindow[dhh_id]->Fill(it.getFrameNr() * 192 - startGate);
163  }
164 
165  if (hrawPxdAdcMapAll) {
166  for (auto& it : m_storeRawAdcs) {
167  int dhh_id;
168  // calculate DHH id from Vxd Id
169  unsigned int layer, ladder, sensor;//, segment;
170  VxdID currentVxdId;
171  currentVxdId = it.getSensorID();
172  layer = currentVxdId.getLayerNumber();
173  ladder = currentVxdId.getLadderNumber();
174  sensor = currentVxdId.getSensorNumber();
175  // segment = currentVxdId.getSegmentNumber();// Frame nr? ... ignore
176  dhh_id = ((layer - 1) << 5) | ((ladder) << 1) | (sensor - 1);
177  if (dhh_id <= 0 || dhh_id >= 64) {
178  B2ERROR("SensorId (DHH ID) out of range: " << dhh_id);
179  continue;
180  }
181 
182  unsigned int chip_offset;
183  chip_offset = it.getChip() * 64;
184  const unsigned char* data = &it.getData()[0];
185  for (int row = 0; row < 768; row++) {
186  for (int col = 0; col < 64; col++) {
187  hrawPxdAdcMapAll->Fill(col + chip_offset + ladder * 300 - 200, 100 + row + 850 * (layer + layer + sensor - 3), *(data++));
188  }
189  }
190  }
191  }
192 }
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 PXD DAQ DHE Status class.
unsigned short getTriggerGate(void) const
get Trigger Gate
unsigned short getDHEID(void) const
Get DHE ID of sensor.
TH2F * hrawPxdAdcMapAll
Histogram Adc 2d hitmap (full frames only) (all pxd)
void initialize() override final
Initialize.
std::string m_storeRawAdcsName
RawAdcs StoreArray name.
TH1F * hrawPxdGateTimeWindow[64]
Histogram raw pixel trigger gate window.
TH1F * hrawPxdPacketSize
Histogram raw packet size.
StoreObjPtr< PXDDAQStatus > m_storeDAQEvtStats
Input array for DAQ Status.
StoreArray< RawPXD > m_storeRawPxdrarray
Storearray for raw data packets
TH1F * hrawPxdHitsCharge[64]
Histogram raw pixel charge.
void defineHisto() override final
Define histograms.
std::string m_storeRawHitsName
PXDRawHits StoreArray name.
std::string m_storeRawPxdrarrayName
RawPXD StoreArray name.
void event() override final
Event.
TH1F * hrawPxdHitsCount
Histogram pixelcount/???
TH2F * hrawPxdHitMapAll
Histogram 2d hitmap (all)
StoreArray< PXDRawAdc > m_storeRawAdcs
Storearray for ADC from full frames
std::string m_histogramDirectoryName
Name of the histogram directory in ROOT file.
TH1F * hrawPxdPackets
Histogram number of raw packets.
TH2F * hrawPxdHitMap[64]
Histogram 2d hitmap.
void beginRun() override final
Begin run.
TH1F * hrawPxdHitTimeWindow[64]
Histogram raw pixel hit "time" window.
StoreArray< PXDRawHit > m_storeRawHits
Storearray for raw pixels
TH2F * hrawPxdChargeMap[64]
Histogram 2d chargemap.
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
bool validSensorID(Belle2::VxdID id) const
Check that id is a valid sensor number.
Definition: GeoCache.cc:52
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
baseType getSensorNumber() const
Get the sensor id.
Definition: VxdID.h:100
baseType getLadderNumber() const
Get the ladder id.
Definition: VxdID.h:98
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
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.