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