Belle II Software  release-06-01-15
PXDDQMCorrModule.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/PXDDQMCorrModule.h"
10 
11 #include <framework/core/HistoModule.h>
12 
13 #include <vector>
14 #include <boost/format.hpp>
15 
16 #include "TH1F.h"
17 #include "TH2F.h"
18 #include "TDirectory.h"
19 
20 using namespace std;
21 using boost::format;
22 using namespace Belle2;
23 
24 //-----------------------------------------------------------------
25 // Register the Module
26 //-----------------------------------------------------------------
27 REG_MODULE(PXDDQMCorr)
28 
29 
30 //-----------------------------------------------------------------
31 // Implementation
32 //-----------------------------------------------------------------
33 
35 {
36  //Set module properties
37  setDescription("PXD DQM Correlation module");
38  setPropertyFlags(c_ParallelProcessingCertified); // specify this flag if you need parallel processing
39  addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms will be placed",
40  std::string("pxd"));
41 }
42 
43 
44 //------------------------------------------------------------------
45 // Function to define histograms
46 //-----------------------------------------------------------------
47 
48 void PXDDQMCorrModule::defineHisto()
49 {
50  // Create a separate histogram directory and cd into it.
51  TDirectory* oldDir = gDirectory;
52  if (m_histogramDirectoryName != "") {
53  oldDir->mkdir(m_histogramDirectoryName.c_str());// do not use return value with ->cd(), its ZERO if dir already exists
54  oldDir->cd(m_histogramDirectoryName.c_str());
55  }
56  //----------------------------------------------------------------
57 
58  int nPixelsU1 = 1;//getInfo(0).getUCells();// we need an actual sensor ID to use that
59  int nPixelsV1 = 3;//getInfo(0).getVCells();
60  int nPixelsU2 = 1;//getInfo(1).getUCells();
61  int nPixelsV2 = 3;//getInfo(1).getVCells();
62  m_CorrelationU = new TH2F("CorrelationU", "Correlation of U;U1/cm;U2/cm", 25, -nPixelsU1, nPixelsU1, 25, nPixelsU2, nPixelsU2);
63  m_CorrelationV = new TH2F("CorrelationV", "Correlation of V;V1/cm;V2/cm", 50, -nPixelsV1, nPixelsV1, 50, nPixelsV2, nPixelsV2);
64  m_DeltaU = new TH1F("DeltaU", "Correlation of U2-U1;Udiff/cm", 100, -nPixelsU1, nPixelsU2);
65  m_DeltaV = new TH1F("DeltaV", "Correlation of V2-V1;Vdiff/cm", 200, -nPixelsV1, nPixelsV2);
66 
67  // cd back to root directory
68  oldDir->cd();
69 }
70 
71 
72 void PXDDQMCorrModule::initialize()
73 {
74  // Register histograms (calls back defineHisto)
75  REG_HISTOGRAM
76 
77  //Register collections
78  m_storeClusters.isRequired(m_storeClustersName);
79 }
80 
81 void PXDDQMCorrModule::beginRun()
82 {
83  // Just to make sure, reset all the histograms.
84  m_CorrelationU->Reset();
85  m_CorrelationV->Reset();
86  m_DeltaU->Reset();
87  m_DeltaV->Reset();
88 }
89 
90 
91 void PXDDQMCorrModule::event()
92 {
93 
94  // If there are no clusters, leave
95  if (!m_storeClusters || !m_storeClusters.getEntries()) return;
96 
97  for (const PXDCluster& cluster1 : m_storeClusters) {
98  int iPlane1 = cluster1.getSensorID().getLayerNumber();
99  if (iPlane1 == 0) {
100  for (const PXDCluster& cluster2 : m_storeClusters) {
101  int iPlane2 = cluster2.getSensorID().getLayerNumber();
102  if (iPlane2 == 1) {
103  m_CorrelationU->Fill(cluster1.getU(), cluster2.getU());
104  m_CorrelationV->Fill(cluster1.getV(), cluster2.getV());
105  m_DeltaU->Fill(cluster2.getU() - cluster1.getU());
106  m_DeltaV->Fill(cluster2.getV() - cluster1.getV());
107  }
108  }
109  }
110  }
111 }
112 
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
The PXD Cluster class This class stores all information about reconstructed PXD clusters The position...
Definition: PXDCluster.h:30
PXD DQM Corr Module.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.