Belle II Software  release-06-01-15
SVDCrossTalkCalibrationsCollectorModule.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 <svd/modules/svdCrossTalkCalibrationsCollector/SVDCrossTalkCalibrationsCollectorModule.h>
10 
11 #include <svd/modules/svdCrossTalkFinder/SVDCrossTalkFinderHelperFunctions.h>
12 
13 #include <vxd/geometry/GeoCache.h>
14 
15 #include <TH1F.h>
16 
17 #include <iostream>
18 
19 
20 using namespace std;
21 using namespace Belle2;
22 
23 
24 REG_MODULE(SVDCrossTalkCalibrationsCollector)
25 
27 {
28  setDescription("This module collects the list of channels exhibiting crossTalk readout.");
29 
30  setPropertyFlags(c_ParallelProcessingCertified);
31 
32  addParam("SVDShaperDigits", m_svdShaperDigitsName,
33  "SVDShaperDigit collection name", string(""));
34 
35  addParam("uSideOccupancyFactor", m_uSideOccupancyFactor,
36  "Multiple of the average occupancy for high occupancy strip classification", 2);
37 
38  addParam("vSideOccupancyFactor", m_vSideOccupancyFactor,
39  "Multiple of the average occupancy for high occupancy strip classification", 2);
40 
41  addParam("nAPVFactor", m_nAPVFactor,
42  "Number of APV chips with at least one high occupancy strips in event required for cross talk flag", 30);
43 
44 }
45 
46 void SVDCrossTalkCalibrationsCollectorModule::prepare()
47 {
48 
49  m_svdShaperDigits.isRequired(m_svdShaperDigitsName);
50 
51 
52  m_histogramTree = new TTree("tree", "tree");
53  m_histogramTree->Branch("hist", "TH1F", &m_hist, 32000, 0);
54  m_histogramTree->Branch("layer", &m_layer, "layer/I");
55  m_histogramTree->Branch("ladder", &m_ladder, "ladder/I");
56  m_histogramTree->Branch("sensor", &m_sensor, "sensor/I");
57  m_histogramTree->Branch("side", &m_side, "side/I");
58 
59  registerObject<TTree>("HTreeCrossTalkCalib", m_histogramTree);
60 
61 
62 }
63 
64 void SVDCrossTalkCalibrationsCollectorModule::startRun()
65 {
66 
67  std::string sensorName;
68  TH1F* sensorHist;
69  //Prepare histograms for payload.
70  VXD::GeoCache& geo = VXD::GeoCache::getInstance();
71  for (auto& layers : geo.getLayers(VXD::SensorInfoBase::SVD)) {
72  for (auto& ladders : geo.getLadders(layers)) {
73  for (auto& sensors : geo.getSensors(ladders)) {
74  for (int side = 0; side <= 1; side++) {
75  occupancyPDFName(sensors, side, sensorName);
76  if (m_sensorHistograms.count(sensorName) == 0) {
77  if (layers.getLayerNumber() == 3 or side == 1) {
78  sensorHist = new TH1F(sensorName.c_str(), "", 768, 0, 768);
79  } else {
80  sensorHist = new TH1F(sensorName.c_str(), "", 512, 0, 512);
81  }
82  m_sensorHistograms[sensorName] = sensorHist;
83 
84  }
85 
86  }
87  }
88  }
89  }
90 
91 }
92 
93 void SVDCrossTalkCalibrationsCollectorModule::collect()
94 {
95  //Vectors to hold event info.
96  vector<std::string> highOccChips_uSide ;
97  vector<std::string> highOccChips_vSide ;
98  vector<std::string> highOccChipsStripNum_uSide ;
99  vector<std::string> highOccChipsStripNum_vSide ;
100  vector<std::string> clusterStrips_uSide ;
101  vector<std::string> clusterStrips_vSide ;
102  vector<std::string> clusterChips_uSide ;
103  vector<int> strips_uSide;
104  vector<int> strips_vSide;
105 
106  //loop over ShaperDigits
107  for (auto& svdShaperDigit : m_svdShaperDigits) {
108  //Remove L3 and +fw sensors, not affected by cross-talk
109  if (svdShaperDigit.getSensorID().getLayerNumber() == 3 or svdShaperDigit.getSensorID().getSensorNumber() == 1) {
110  continue;
111  }
112 
113  int side = svdShaperDigit.isUStrip();
114  std::string sensorName;
115  occupancyPDFName(svdShaperDigit.getSensorID(), side, sensorName);
116 
117  double sensorAverage = 0.;
118  calculateAverage(svdShaperDigit.getSensorID(), sensorAverage, side);
119  int stripID = svdShaperDigit.getCellID();
120  std::string sensorStripNum = sensorName + "." + std::to_string(stripID);
121  double stripOccupancy = m_OccupancyCal.getOccupancy(svdShaperDigit.getSensorID(), side, stripID);
122  //Clustering only works assuming digits are ordered.//
123  if (side == 1 && stripOccupancy > (m_uSideOccupancyFactor * sensorAverage)) {
124 
125  int adjacentStrip = 0;
126  if (!strips_uSide.empty()) {
127  adjacentStrip = stripID - strips_uSide.back();
128  }
129  strips_uSide.push_back(stripID);
130  if (!highOccChips_uSide.empty() && sensorName == highOccChips_uSide.back() && adjacentStrip == 1) {
131  clusterChips_uSide.push_back(sensorName);
132  if (clusterStrips_uSide.empty() || clusterStrips_uSide.back() != sensorStripNum) {
133  clusterStrips_uSide.push_back(highOccChipsStripNum_uSide.back());
134  }
135  clusterStrips_uSide.push_back(sensorStripNum);
136  }
137  highOccChipsStripNum_uSide.push_back(sensorStripNum);
138  highOccChips_uSide.push_back(sensorName);
139  }
140 
141  if (side == 0 && stripOccupancy > (m_vSideOccupancyFactor * sensorAverage)) {
142  int adjacentStrip = 0;
143  if (!strips_vSide.empty()) {
144  adjacentStrip = stripID - strips_vSide.back();
145  }
146  strips_vSide.push_back(stripID);
147  if (!highOccChips_vSide.empty() && sensorName == highOccChips_vSide.back() && adjacentStrip == 1) {
148  if (clusterStrips_vSide.empty() || clusterStrips_vSide.back() != sensorStripNum) {
149  clusterStrips_vSide.push_back(highOccChipsStripNum_vSide.back());
150  }
151  clusterStrips_vSide.push_back(sensorStripNum);
152  }
153  highOccChipsStripNum_vSide.push_back(sensorStripNum);
154  highOccChips_vSide.push_back(sensorName);
155  }
156 
157  } //ShaperDigit loop
158 
159  std::sort(clusterChips_uSide.begin(), clusterChips_uSide.end());
160  clusterChips_uSide.erase(unique(clusterChips_uSide.begin(), clusterChips_uSide.end()), clusterChips_uSide.end());
161  int numberOfClusterChips = std::distance(clusterChips_uSide.begin(), std::unique(clusterChips_uSide.begin(),
162  clusterChips_uSide.end()));
163 
164 // Crosstalk events flagged using u-side
165  if (numberOfClusterChips > m_nAPVFactor) {
166  for (auto& svdShaperDigit : m_svdShaperDigits) {
167  std::string sensorID = svdShaperDigit.getSensorID();
168  std::string digitID = sensorID + "." + std::to_string(svdShaperDigit.isUStrip());
169  std::string stripID = digitID + "." + std::to_string(svdShaperDigit.getCellID());
170  if (std::find(clusterStrips_uSide.begin(), clusterStrips_uSide.end(), stripID) != clusterStrips_uSide.end()) {
171  std::string sensorName;
172  occupancyPDFName(svdShaperDigit.getSensorID(), svdShaperDigit.isUStrip(), sensorName);
173  auto xTalkStrip = m_sensorHistograms.at(sensorName);
174  xTalkStrip->Fill(svdShaperDigit.getCellID(), 1);
175  }
176  if (std::find(clusterStrips_vSide.begin(), clusterStrips_vSide.end(), stripID) != clusterStrips_vSide.end()) {
177  std::string sensorName;
178  occupancyPDFName(svdShaperDigit.getSensorID(), svdShaperDigit.isUStrip(), sensorName);
179  auto xTalkStrip = m_sensorHistograms.at(sensorName);
180  xTalkStrip->Fill(svdShaperDigit.getCellID(), 1);
181  }
182 
183 
184  }//shaper digit loop
185  }
186 } //Collector loop
187 
188 
189 void SVDCrossTalkCalibrationsCollectorModule::finish()
190 {}
191 
192 void SVDCrossTalkCalibrationsCollectorModule::closeRun()
193 {
194  std::string sensorName;
195  VXD::GeoCache& geo = VXD::GeoCache::getInstance();
196  for (auto& layers : geo.getLayers(VXD::SensorInfoBase::SVD)) {
197  for (auto& ladders : geo.getLadders(layers)) {
198  for (auto& sensors : geo.getSensors(ladders)) {
199  for (int side = 0; side <= 1; side++) {
200 
201  occupancyPDFName(sensors, side, sensorName);
202  auto sensorOnMap = m_sensorHistograms.at(sensorName);
203  m_hist = sensorOnMap;
204  m_layer = sensors.getLayerNumber();
205  m_ladder = sensors.getLadderNumber();
206  m_sensor = sensors.getSensorNumber();
207  m_side = side;
208 
209  getObjectPtr<TTree>("HTreeCrossTalkCalib")->Fill();
210 
211  sensorOnMap->Delete();
212  }
213  }
214  }
215  }
216 
217 }
218 
219 
220 void SVDCrossTalkCalibrationsCollectorModule::calculateAverage(const VxdID& sensorID, double& mean, int side)
221 {
222  double nBins = 0;
223  if (side == 1) nBins = 768; //U-side 768 channels
224  else nBins = 512;
225  double count = 0;
226  for (int i = 0; i < nBins; i++) {
227  count += m_OccupancyCal.getOccupancy(sensorID, side, i);
228  }
229  mean = count / nBins;
230 }
231 
Calibration collector module base class.
Collector module used to create the histograms needed for the cross talk calibration.
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:39
const std::set< Belle2::VxdID > getLayers(SensorInfoBase::SensorType sensortype=SensorInfoBase::VXD)
Return a set of all known Layers.
Definition: GeoCache.cc:175
const std::set< Belle2::VxdID > & getSensors(Belle2::VxdID ladder) const
Return a set of all sensor IDs belonging to a given ladder.
Definition: GeoCache.cc:203
const std::set< Belle2::VxdID > & getLadders(Belle2::VxdID layer) const
Return a set of all ladder IDs belonging to a given layer.
Definition: GeoCache.cc:192
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
void occupancyPDFName(const VxdID &sensor, int side, std::string &PDFName)
Function to maintain common naming convention between calibration occupancy file generation and occup...
Abstract base class for different kinds of events.