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