Belle II Software  release-06-02-00
PXDROIPlotModule.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/pxdHelper/PXDROIPlotModule.h>
10 #include <TBox.h>
11 #include <TCanvas.h>
12 #include <TStyle.h>
13 #include <TH2.h>
14 #include <boost/format.hpp>
15 
16 using namespace std;
17 using namespace Belle2;
18 using namespace Belle2::PXD;
19 
20 //-----------------------------------------------------------------
21 // Register the Module
22 //-----------------------------------------------------------------
23 REG_MODULE(PXDROIPlot)
24 
25 //-----------------------------------------------------------------
26 // Implementation
27 //-----------------------------------------------------------------
28 
30 {
31  //Set module properties
32  setDescription("Plot ROIs on PXD Hit/Charge Maps and write pictures");
33  setPropertyFlags(c_ParallelProcessingCertified);
34  addParam("PXDRawHitsName", m_storeRawHitsName, "The name of the StoreArray of PXDRawHits to be processed", string(""));
35  // we do not want to have a default name, as the user should be very sure what he is plotting!
36  addParam("ROIsName", m_ROIsName, "name of the list of ROIs (plotted in red)", string("__unusedROIs__"));
37  addParam("DCROIsName", m_DCROIsName, "name of the list of DatCon ROIs (optional, plotted in blue)", string("__unusedROIs__"));
38  addParam("HLTROIsName", m_HLTROIsName, "name of the list of HLT ROIs (optional, plotted in green)", string("__unusedROIs__"));
39  addParam("Prefix", m_prefix, "prefix for plots", string(""));
40 }
41 
42 void PXDROIPlotModule::initialize()
43 {
44  m_eventMetaData.isRequired();
45  m_storeRawHits.isRequired(m_storeRawHitsName);
46  // all ROIs are optional, allowing to use always the same color for comparison even if one cathegory is not there
47  m_storeROIs.isOptional(m_ROIsName);
48  m_storeDCROIs.isOptional(m_DCROIsName);
49  m_storeHLTROIs.isOptional(m_HLTROIsName);
50 }
51 
52 void PXDROIPlotModule::event()
53 {
54  unsigned int evtNr = m_eventMetaData->getEvent();
55  unsigned int evtRun = m_eventMetaData->getRun();
56  map <VxdID, bool> flag; // flag sensors with ROIs
57  map <VxdID, vector <ROIid>> list;
58  map <VxdID, vector <ROIid>> listDC;
59  map <VxdID, vector <ROIid>> listHLT;
60  string histoname;
61  TCanvas* c1;
62  gStyle->SetPalette(55);
63  gStyle->SetOptStat(0);
64  c1 = new TCanvas("c1", "c1", 1200, 600);
65  TH2F* h = new TH2F("ChargemapWithROIs", "PXD Module Chargemap;VCell;UCell", 768, 0, 768, 250, 0, 250);
66  h->SetContour(100);
67  for (auto& it : m_storeROIs) {
68  list[it.getSensorID()].push_back(it);
69  flag[it.getSensorID()] = true;
70  }
71  for (auto& it : m_storeDCROIs) {
72  listDC[it.getSensorID()].push_back(it);
73  flag[it.getSensorID()] = true;
74  }
75  for (auto& it : m_storeHLTROIs) {
76  listHLT[it.getSensorID()].push_back(it);
77  flag[it.getSensorID()] = true;
78  }
79 
80  for (auto& f : flag) {
81  c1->Clear();
82  c1->cd();
83  h->Reset();
84 
85  VxdID currentVxdId = f.first;
86  histoname = m_prefix + boost::str(boost::format("Run_%d_Evt_%d_") % evtRun % evtNr) + string(currentVxdId);
87  h->SetTitle(histoname.data());
88 
89  for (auto& pix : m_storeRawHits) {
90 
91  if (currentVxdId != pix.getSensorID()) continue;
92  h->Fill(pix.getRow(), pix.getColumn(), pix.getCharge());
93 
94  }
95 
96  h->Draw("colz");
97 
98  // Alpha Blending only works in ROOT's Batch Mode, but LineStyle won't.
99  // LineStyle works in normal mode, but Alpha won't
100  for (auto& it : list[currentVxdId]) {
101  TBox* b;
102  b = new TBox(it.getMinVid(), it.getMinUid(), it.getMaxVid(), it.getMaxUid());
103  b->SetLineColorAlpha(kRed, 0.3);
104  b->SetLineWidth(3);
105  b->SetFillStyle(0);// Hollow
106  b->Draw();
107  }
108  // we move the other box by more than half a bin. this is needed as alpha seems not to work
109  for (auto& it : listDC[currentVxdId]) {
110  TBox* b;
111  b = new TBox(it.getMinVid() + 0.7, it.getMinUid() + 0.7, it.getMaxVid() + 0.7, it.getMaxUid() + 0.7);
112  b->SetLineColorAlpha(kBlue, 0.3);
113  b->SetLineWidth(2);
114  b->SetLineStyle(2);
115  b->SetFillStyle(0);// Hollow
116  b->Draw();
117  }
118  // we move the other box by half a bin. this is needed as alpha seems not to work, in addition we use a dashed style
119  // dashed style doesnt work with png export, thus if all ROIs are identical, lines might overlap completely
120  for (auto& it : listHLT[currentVxdId]) {
121  TBox* b;
122  b = new TBox(it.getMinVid() - 0.7, it.getMinUid() - 0.7, it.getMaxVid() - 0.7, it.getMaxUid() - 0.7);
123  b->SetLineColorAlpha(kGreen, 0.3);
124  b->SetLineWidth(2);
125  b->SetLineStyle(7);
126  b->SetFillStyle(0);// Hollow
127  b->Draw();
128  }
129  c1->Print((histoname + ".png").data());
130  c1->Print((histoname + ".root").data());
131  }
132  delete h;
133  delete c1;
134 }
Base class for Modules.
Definition: Module.h:72
Plot each event with ROI and Pixels.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
#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.
Abstract base class for different kinds of events.