Belle II Software  release-08-01-10
PXDRawDQMChipsModule.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/PXDRawDQMChipsModule.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(PXDRawDQMChips);
27 
28 //-----------------------------------------------------------------
29 // Implementation
30 //-----------------------------------------------------------------
31 
32 PXDRawDQMChipsModule::PXDRawDQMChipsModule() : HistoModule(), m_storeRawHits()
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("PXDRawHitsName", m_storeRawHitsName, "The name of the StoreArray of PXDRawHits to be processed", string(""));
41 }
42 
44 {
45  // Create a separate histogram directory and cd into it.
46  TDirectory* oldDir = gDirectory;
47  if (m_histogramDirectoryName != "") {
48  oldDir->mkdir(m_histogramDirectoryName.c_str());
49  oldDir->cd(m_histogramDirectoryName.c_str());
50  }
51 
52 
53  for (auto i = 0; i < 64; i++) {
54  auto layer = (((i >> 5) & 0x1) + 1);
55  auto ladder = ((i >> 1) & 0xF);
56  auto sensor = ((i & 0x1) + 1);
57 
58  // Check if sensor exist
59  if (Belle2::VXD::GeoCache::getInstance().validSensorID(Belle2::VxdID(layer, ladder, sensor))) {
60  for (auto j = 0; j < eNumSwitcher; j++) {
61  for (auto k = 0; k < eNumDCD; k++) {
62  string s = str(format("Sensor %d:%d:%d (DHH ID %02Xh) Switcher %d DCD %d") % layer % ladder % sensor % i % j % k);
63  string s2 = str(format("_%d.%d.%d_%d_%d") % layer % ladder % sensor % j % k);
64 
65  hrawPxdHitsCount[i][j][k] = new TH1F(("hrawPxdCount" + s2).c_str(), ("Pxd Raw Count " + s + ";Nr per Event").c_str(), 8192, 0,
66  8192);
67  hrawPxdHitsCharge[i][j][k] = new TH1F(("hrawPxdHitsCharge" + s2).c_str(),
68  ("Pxd Raw Hit Charge, " + s + ";Charge").c_str(), 256, 0, 256);
69  }
70  }
71  } else {
72  for (auto j = 0; j < eNumSwitcher; j++) {
73  for (auto k = 0; k < eNumDCD; k++) {
74  hrawPxdHitsCount[i][j][k] = NULL;
75  hrawPxdHitsCharge[i][j][k] = NULL;
76  }
77  }
78  }
79  }
80 
81  // cd back to root directory
82  oldDir->cd();
83 }
84 
86 {
87  REG_HISTOGRAM
89 }
90 
92 {
93  // Just to make sure, reset all the histograms.
94  for (auto i = 0; i < eNumSensors; i++) {
95  for (auto j = 0; j < eNumSwitcher; j++) {
96  for (auto k = 0; k < eNumDCD; k++) {
97  if (hrawPxdHitsCount[i][j][k]) hrawPxdHitsCount[i][j][k]->Reset();
98  if (hrawPxdHitsCharge[i][j][k]) hrawPxdHitsCharge[i][j][k]->Reset();
99  }
100  }
101  }
102 }
103 
105 {
106  uint nhits[eNumSensors][eNumSwitcher][eNumDCD] = {};
107 
108  for (auto& it : m_storeRawHits) {
109  int dhh_id;
110  // calculate DHH id from Vxd Id
111  unsigned int layer, ladder, sensor;//, segment;
112  VxdID currentVxdId;
113  currentVxdId = it.getSensorID();
114  layer = currentVxdId.getLayerNumber();
115  ladder = currentVxdId.getLadderNumber();
116  sensor = currentVxdId.getSensorNumber();
117  // segment = currentVxdId.getSegmentNumber();// Frame nr? ... ignore
118  dhh_id = ((layer - 1) << 5) | ((ladder) << 1) | (sensor - 1);
119  if (dhh_id <= 0 || dhh_id >= 64) {
120  B2ERROR("SensorId (DHH ID) out of range: " << dhh_id);
121  continue;
122  }
123  int switcher = it.getRow() / 128; // TODO here must be: function for REAL inverse mapping
124  int dcd = it.getColumn() / 64; // TODO here must be: function for REAL inverse mapping
125  // TODO check switcher 0-6? DCD 0-4?
126  nhits[dhh_id][switcher][dcd]++;
127  if (hrawPxdHitsCharge[dhh_id][switcher][dcd]) hrawPxdHitsCharge[dhh_id][switcher][dcd]->Fill(it.getCharge());
128  }
129  for (auto i = 0; i < eNumSensors; i++) {
130  for (auto j = 0; j < eNumSwitcher; j++) {
131  for (auto k = 0; k < eNumDCD; k++) {
132  if (hrawPxdHitsCount[i][j][k]) hrawPxdHitsCount[i][j][k]->Fill(nhits[i][j][k]);
133  }
134  }
135  }
136 }
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
virtual void initialize() override
Initialize.
virtual void event() override
Event.
TH1F * hrawPxdHitsCharge[eNumSensors][eNumSwitcher][eNumDCD]
Histogram raw pixel charge.
std::string m_storeRawHitsName
PXDRawHits StoreArray name.
TH1F * hrawPxdHitsCount[eNumSensors][eNumSwitcher][eNumDCD]
Histogram pixelcount/???
virtual void beginRun() override
Begin run.
std::string m_histogramDirectoryName
Name of the histogram directory in ROOT file.
StoreArray< PXDRawHit > m_storeRawHits
Storearray for raw pixels
virtual void defineHisto() override
Define histograms.
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.