Belle II Software  release-08-01-10
SVDROIDQMModule.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 <tracking/modules/svdROIFinder/SVDROIDQMModule.h>
10 
11 #include <TDirectory.h>
12 #include <TH2F.h>
13 
14 using namespace Belle2;
15 
16 //-----------------------------------------------------------------
17 // Register the Module
18 //-----------------------------------------------------------------
19 REG_MODULE(SVDROIDQM);
20 
21 //-----------------------------------------------------------------
22 // Implementation
23 //-----------------------------------------------------------------
24 
26  : HistoModule()
27  , m_InterDir(nullptr)
28  , m_ROIDir(nullptr)
29  , hInterDictionary(172, [](const Belle2::VxdID & vxdid) {return (size_t)vxdid.getID(); })
30 , hROIDictionary(172, [](const Belle2::VxdID& vxdid) {return (size_t)vxdid.getID(); })
31 , hROIDictionaryEvt(172, [](const Belle2::VxdID& vxdid) {return (size_t)vxdid.getID(); })
32 , m_numModules(0)
33 , hnROIs(nullptr)
34 , hnInter(nullptr)
35 , harea(nullptr)
36 , hredFactor(nullptr)
37 , hCellU(nullptr)
38 , hCellV(nullptr)
39 , n_events(0)
40 {
41  //Set module properties
42  setDescription("Monitor of the ROIs creation on HLT");
43  setPropertyFlags(c_ParallelProcessingCertified);
44 
45  addParam("SVDShaperDigitsName", m_SVDShaperDigitsName,
46  "name of the list of SVDShaperDigits", std::string(""));
47  addParam("SVDRecoDigitsName", m_SVDRecoDigitsName,
48  "name of the list of SVDRecoDigits", std::string(""));
49  addParam("SVDClustersName", m_SVDClustersName,
50  "name of the list of SVDClusters", std::string(""));
51 
52  addParam("InterceptsName", m_InterceptsName,
53  "name of the list of interceptions", std::string(""));
54 
55  addParam("ROIsName", m_ROIsName,
56  "name of the list of ROIs", std::string(""));
57 
58  addParam("specificLayer", m_specificLayer,
59  "Layer number, if you want the plots only for a specific SVD layer. If it is not a SVD layer (3, 4, 5, 6) than the plots for all SVD layers are produced. Default is (-1), i.e. plots for all SVD layers are produced.",
60  m_specificLayer);
61 
62  addParam("plotRecoDigits", m_plotRecoDigits,
63  "Set true to produce the plots for RecoDigits (false by default)", m_plotRecoDigits);
64 
65 }
66 
68 {
69 
70  // Create a separate histogram directory and cd into it.
71  TDirectory* oldDir = gDirectory;
72  TDirectory* roiDir = oldDir->mkdir("SVDROIs");
73  m_InterDir = roiDir->mkdir("intercept");
74  m_ROIDir = roiDir->mkdir("roi");
75 
76  hCellU = new TH1F("hCellU", "CellID U", 769, -0.5, 768.5);
77  hCellU->GetXaxis()->SetTitle("U cell ID");
78  hCellV = new TH1F("hCellV", "CellID V", 769, -0.5, 768.5);
79  hCellV->GetXaxis()->SetTitle("V cell ID");
80 
81  m_InterDir->cd();
82  hnInter = new TH1F("hnInter", "number of intercepts", 100, 0, 100);
83 
84  m_ROIDir->cd();
85  hnROIs = new TH1F("hnROIs", "number of ROIs", 100, 0, 100);
86  harea = new TH1F("harea", "ROIs area", 100, 0, 100000);
87  hredFactor = new TH1F("hredFactor", "ROI reduction factor", 1000, 0, 1);
88 
89 
91 
92  oldDir->cd();
93 
94 }
95 
97 {
98  REG_HISTOGRAM
99 
103  m_ROIs.isRequired(m_ROIsName);
104  m_Intercepts.isRequired(m_InterceptsName);
105 
106  n_events = 0;
107 
108 
109 }
110 
112 {
113 
114  n_events++;
115 
116  for (auto& it : m_SVDShaperDigits)
117  if (it.isUStrip())
118  hCellU->Fill(it.getCellID());
119  else
120  hCellV->Fill(it.getCellID());
121 
122  hnInter->Fill(m_Intercepts.getEntries());
123 
124  for (auto& it : m_Intercepts)
126 
127 
128  for (auto it = hROIDictionaryEvt.begin(); it != hROIDictionaryEvt.end(); ++it)
129  (it->second).value = 0;
130 
131  int ROIarea = 0;
132  double redFactor = 0;
133 
134  for (auto& it : m_ROIs) {
135 
136  fillSensorROIHistos(&it);
137 
138  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
139  const int nPixelsU = aSensorInfo.getUCells();
140  const int nPixelsV = aSensorInfo.getVCells();
141 
142  int minU = it.getMinUid();
143  int minV = it.getMinVid();
144  int maxU = it.getMaxUid();
145  int maxV = it.getMaxVid();
146 
147  int tmpROIarea = (maxU - minU) * (maxV - minV);
148  ROIarea += tmpROIarea;
149  redFactor += (double)tmpROIarea / (nPixelsU * nPixelsV * m_numModules);
150 
151  }
152 
153  hnROIs->Fill(m_ROIs.getEntries());
154 
155  harea->Fill((double)ROIarea);
156 
157  hredFactor->Fill((double)redFactor);
158 
159 
160  for (auto it = hROIDictionaryEvt.begin(); it != hROIDictionaryEvt.end(); ++it) {
161  ROIHistoAccumulateAndFill aROIHistoAccumulateAndFill = it->second;
162  aROIHistoAccumulateAndFill.fill(aROIHistoAccumulateAndFill.hPtr, aROIHistoAccumulateAndFill.value);
163  }
164 
165 }
166 
168 {
169 
170  // VXD::GeoCache& aGeometry = VXD::GeoCache::getInstance();
171 
172  std::string name; //name of the histogram
173  std::string title; //title of the histogram
174  TH2F* tmp2D; //temporary 2D histo used to set axis title
175  TH1F* tmp1D; //temporary 1D histo used to set axis title
176 
177  m_numModules = 0;
178 
179  std::set<Belle2::VxdID> svdLayers = m_geoCache.getLayers(VXD::SensorInfoBase::SVD);
180  std::set<Belle2::VxdID>::iterator itSvdLayers = svdLayers.begin();
181 
182  if (m_specificLayer >= 3 && m_specificLayer <= 6) {
183  B2INFO("Producing plots for layer: " << m_specificLayer);
184  svdLayers.clear();
185  svdLayers.insert(Belle2::VxdID(m_specificLayer, 0, 0));
186  itSvdLayers = svdLayers.begin();
187  } else {
188  B2INFO("No specific SVD layer (3,4,5,6) selected (m_specificLayer = " << m_specificLayer <<
189  "). Producing plots for all SVD layers.");
190  }
191 
192  while (itSvdLayers != svdLayers.end()) {
193 
194  std::set<Belle2::VxdID> svdLadders = m_geoCache.getLadders(*itSvdLayers);
195  std::set<Belle2::VxdID>::iterator itSvdLadders = svdLadders.begin();
196 
197  while (itSvdLadders != svdLadders.end()) {
198 
199  std::set<Belle2::VxdID> svdSensors = m_geoCache.getSensors(*itSvdLadders);
200  std::set<Belle2::VxdID>::iterator itSvdSensors = svdSensors.begin();
201 
202  while (itSvdSensors != svdSensors.end()) {
203 
204  m_numModules++; //counting the total number of modules
205 
206  const VXD::SensorInfoBase& wSensorInfo = m_geoCache.getSensorInfo(*itSvdSensors);
207 
208  const int nPixelsU = wSensorInfo.getUCells();
209  const int nPixelsV = wSensorInfo.getVCells();
210  std::string sensorid = std::to_string(itSvdSensors->getLayerNumber()) + "_" + std::to_string(
211  itSvdSensors->getLadderNumber()) + "_" +
212  std::to_string(itSvdSensors->getSensorNumber());
213 
214 
215  // ------ HISTOGRAMS WITH AN ACCUMULATE PER ROI AND A FILL PER EVENT -------
216  m_ROIDir->cd();
217 
218  name = "hNROIs_" + sensorid;
219  title = "number of ROIs for sensor " + sensorid;
220  double value = 0;
222  new TH1F(name.c_str(), title.c_str(), 25, 0, 25),
223  [](const ROIid*, double & val) {val++;},
224  [](TH1 * hPtr, double & val) { hPtr->Fill(val); },
225  value
226  };
227  hROIDictionaryEvt.insert(std::pair< Belle2::VxdID, ROIHistoAccumulateAndFill& > ((Belle2::VxdID)*itSvdSensors, *aHAAF));
228 
229 
230 
231 
232  // ------ HISTOGRAMS WITH A FILL PER INTERCEPT -------
233  m_InterDir->cd();
234 
235  // coor U and V
236  name = "hCoorU_" + sensorid;
237  title = "U coordinate of the extrapolation in U for sensor " + sensorid;
238  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
239  (
240  (Belle2::VxdID)*itSvdSensors,
242  new TH1F(name.c_str(), title.c_str(), 100, -5, 5),
243  [](TH1 * hPtr, const SVDIntercept * inter) { hPtr->Fill(inter->getCoorU()); }
244  )
245  )
246  );
247 
248  name = "hCoorV_" + sensorid;
249  title = "V coordinate of the extrapolation in V for sensor " + sensorid;
250  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
251  (
252  (Belle2::VxdID)*itSvdSensors,
254  new TH1F(name.c_str(), title.c_str(), 100, -5, 5),
255  [](TH1 * hPtr, const SVDIntercept * inter) { hPtr->Fill(inter->getCoorV()); }
256  )
257  )
258  );
259 
260  // Intercept U vs V coordinate
261  name = "hCoorU_vs_CoorV_" + sensorid;
262  title = "U vs V intercept (cm) " + sensorid;
263  tmp2D = new TH2F(name.c_str(), title.c_str(), 100, -5, 5, 100, -5, 5);
264  tmp2D->GetXaxis()->SetTitle("intercept U coor (cm)");
265  tmp2D->GetYaxis()->SetTitle("intercept V coor (cm)");
266  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
267  (
268  (Belle2::VxdID)*itSvdSensors,
270  tmp2D,
271  [](TH1 * hPtr, const SVDIntercept * inter) { hPtr->Fill(inter->getCoorU(), inter->getCoorV()); }
272  )
273  )
274  );
275 
276 
277  // sigma U and V
278  name = "hStatErrU_" + sensorid;
279  title = "stat error of the extrapolation in U for sensor " + sensorid;
280  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
281  (
282  (Belle2::VxdID)*itSvdSensors,
284  new TH1F(name.c_str(), title.c_str(), 100, 0, 0.35),
285  [](TH1 * hPtr, const SVDIntercept * inter) { hPtr->Fill(inter->getSigmaU()); }
286  )
287  )
288  );
289  name = "hStatErrV_" + sensorid;
290  title = "stat error of the extrapolation in V for sensor " + sensorid;
291  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
292  (
293  (Belle2::VxdID)*itSvdSensors,
295  new TH1F(name.c_str(), title.c_str(), 100, 0, 0.35),
296  [](TH1 * hPtr, const SVDIntercept * inter) { hPtr->Fill(inter->getSigmaV()); }
297  )
298  )
299  );
300 
301  //1D residuals
302  name = "hResidU_" + sensorid;
303  title = "U residuals = intercept - digit, for sensor " + sensorid;
304  tmp1D = new TH1F(name.c_str(), title.c_str(), 1000, -5, 5);
305  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
306  (
307  (Belle2::VxdID)*itSvdSensors,
309  tmp1D,
310  [this](TH1 * hPtr, const SVDIntercept * inter) {
311 
312  for (auto& it : this->m_SVDShaperDigits)
313  if ((int)it.getSensorID() == (int)inter->getSensorID()) {
314  if (it.isUStrip()) {
315  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
316  hPtr->Fill(inter->getCoorU() - aSensorInfo.getUCellPosition(it.getCellID()));
317  }
318  }
319  }
320  )
321  )
322  );
323 
324  name = "hResidV_" + sensorid;
325  title = "V residuals = intercept - digit, for sensor " + sensorid;
326  tmp1D = new TH1F(name.c_str(), title.c_str(), 1000, -5, 5);
327  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
328  (
329  (Belle2::VxdID)*itSvdSensors,
331  tmp1D,
332  [this](TH1 * hPtr, const SVDIntercept * inter) {
333 
334  for (auto& it : this->m_SVDShaperDigits)
335  if ((int)it.getSensorID() == (int)inter->getSensorID()) {
336  if (!it.isUStrip()) {
337  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
338  hPtr->Fill(inter->getCoorV() - aSensorInfo.getVCellPosition(it.getCellID()));
339  }
340  }
341  }
342  )
343  )
344  );
345 
346 
347  //residual U,V vs coordinate U,V
348  name = "hResidU_vs_CoorU_" + sensorid;
349  title = "U residual (cm) vs coor U (cm) " + sensorid;
350  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -5, 5, 1000, -5, 5);
351  tmp2D->GetYaxis()->SetTitle("U resid (cm)");
352  tmp2D->GetXaxis()->SetTitle("U coor (cm)");
353  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
354  (
355  (Belle2::VxdID)*itSvdSensors,
357  tmp2D,
358  [this](TH1 * hPtr, const SVDIntercept * inter) {
359 
360  for (auto& it : this->m_SVDShaperDigits)
361  if (((int)it.getSensorID() == (int)inter->getSensorID()) && it.isUStrip()) {
362  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
363  double resid = inter->getCoorU() - aSensorInfo.getUCellPosition(it.getCellID());
364  hPtr->Fill(inter->getCoorU(), resid);
365  }
366  }
367  )
368  )
369  );
370 
371  name = "hResidV_vs_CoorV_" + sensorid;
372  title = "V residual (cm) vs coor V (cm) " + sensorid;
373  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -5, 5, 1000, -5, 5);
374  tmp2D->GetYaxis()->SetTitle("V resid (cm)");
375  tmp2D->GetXaxis()->SetTitle("V coor (cm)");
376  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
377  (
378  (Belle2::VxdID)*itSvdSensors,
380  tmp2D,
381  [this](TH1 * hPtr, const SVDIntercept * inter) {
382 
383  for (auto& it : this->m_SVDShaperDigits)
384  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (!it.isUStrip())) {
385  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
386  double resid = inter->getCoorV() - aSensorInfo.getVCellPosition(it.getCellID());
387  hPtr->Fill(inter->getCoorV(), resid);
388  }
389  }
390  )
391  )
392  );
393 
394 
395  // RecoDigits
396  //residual vs charge
397  if (m_plotRecoDigits) {
398  name = "hResidU_vs_charge_" + sensorid;
399  title = "U residual (cm) vs charge " + sensorid;
400  tmp2D = new TH2F(name.c_str(), title.c_str(), 250, 0, 250, 100, -5, 5);
401  tmp2D->GetYaxis()->SetTitle("U resid (cm)");
402  tmp2D->GetXaxis()->SetTitle("charge");
403  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
404  (
405  (Belle2::VxdID)*itSvdSensors,
407  tmp2D,
408  [this](TH1 * hPtr, const SVDIntercept * inter) {
409 
410  for (auto& it : this->m_SVDRecoDigits)
411  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (it.isUStrip())) {
412  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
413  double resid = inter->getCoorU() - aSensorInfo.getUCellPosition(it.getCellID());
414  hPtr->Fill(it.getCharge(), resid);
415  }
416  }
417  )
418  )
419  );
420 
421  name = "hResidV_vs_charge_" + sensorid;
422  title = "V residual (cm) vs charge " + sensorid;
423  tmp2D = new TH2F(name.c_str(), title.c_str(), 250, 0, 250, 100, -5, 5);
424  tmp2D->GetYaxis()->SetTitle("V resid (cm)");
425  tmp2D->GetXaxis()->SetTitle("charge");
426  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
427  (
428  (Belle2::VxdID)*itSvdSensors,
430  tmp2D,
431  [this](TH1 * hPtr, const SVDIntercept * inter) {
432 
433  for (auto& it : this->m_SVDRecoDigits)
434  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (!it.isUStrip())) {
435  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
436  double resid = inter->getCoorV() - aSensorInfo.getVCellPosition(it.getCellID());
437  hPtr->Fill(it.getCharge(), resid);
438  }
439  }
440  )
441  )
442  );
443  }
444 
445  // 1D residual for clusters
446  name = "hClusterResidU_" + sensorid;
447  title = "Cluster U residuals = intercept - cluster, for sensor " + sensorid;
448  tmp1D = new TH1F(name.c_str(), title.c_str(), 1000, -5, 5);
449  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
450  (
451  (Belle2::VxdID)*itSvdSensors,
453  tmp1D,
454  [this](TH1 * hPtr, const SVDIntercept * inter) {
455 
456  for (auto& it : this->m_SVDClusters)
457  if ((int)it.getSensorID() == (int)inter->getSensorID()) {
458  if (it.isUCluster()) {
459  hPtr->Fill(inter->getCoorU() - it.getPosition(inter->getCoorV()));
460  }
461  }
462  }
463  )
464  )
465  );
466 
467  name = "hClusterResidV_" + sensorid;
468  title = "Cluster V residuals = intercept - cluster, for sensor " + sensorid;
469  tmp1D = new TH1F(name.c_str(), title.c_str(), 1000, -5, 5);
470  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
471  (
472  (Belle2::VxdID)*itSvdSensors,
474  tmp1D,
475  [this](TH1 * hPtr, const SVDIntercept * inter) {
476 
477  for (auto& it : this->m_SVDClusters)
478  if ((int)it.getSensorID() == (int)inter->getSensorID()) {
479  if (!it.isUCluster()) {
480  hPtr->Fill(inter->getCoorV() - it.getPosition());
481  }
482  }
483  }
484  )
485  )
486  );
487 
488  //residual U,V vs coordinate U,V for clusters
489  name = "hClusterResidU_vs_CoorU_" + sensorid;
490  title = "Cluster U residual (cm) vs coor U (cm) " + sensorid;
491  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -5, 5, 1000, -5, 5);
492  tmp2D->GetYaxis()->SetTitle("Cluster U resid (cm)");
493  tmp2D->GetXaxis()->SetTitle("U coor (cm)");
494  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
495  (
496  (Belle2::VxdID)*itSvdSensors,
497  InterHistoAndFill(
498  tmp2D,
499  [this](TH1 * hPtr, const SVDIntercept * inter) {
500 
501  for (auto& it : this->m_SVDClusters)
502  if (((int)it.getSensorID() == (int)inter->getSensorID()) && it.isUCluster()) {
503  double resid = inter->getCoorU() - it.getPosition(inter->getCoorV());
504  hPtr->Fill(inter->getCoorU(), resid);
505  }
506  }
507  )
508  )
509  );
510 
511  name = "hClusterResidV_vs_CoorV_" + sensorid;
512  title = "Cluster V residual (cm) vs coor V (cm) " + sensorid;
513  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -5, 5, 1000, -5, 5);
514  tmp2D->GetYaxis()->SetTitle("Cluster V resid (cm)");
515  tmp2D->GetXaxis()->SetTitle("V coor (cm)");
516  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
517  (
518  (Belle2::VxdID)*itSvdSensors,
519  InterHistoAndFill(
520  tmp2D,
521  [this](TH1 * hPtr, const SVDIntercept * inter) {
522 
523  for (auto& it : this->m_SVDClusters)
524  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (!it.isUCluster())) {
525  double resid = inter->getCoorV() - it.getPosition();
526  hPtr->Fill(inter->getCoorV(), resid);
527  }
528  }
529  )
530  )
531  );
532 
533 
534  //residual vs charge for clusters
535  name = "hClusterResidU_vs_charge_" + sensorid;
536  title = "Cluster U residual (cm) vs charge " + sensorid;
537  tmp2D = new TH2F(name.c_str(), title.c_str(), 250, 0, 250, 100, -5, 5);
538  tmp2D->GetYaxis()->SetTitle("U resid (cm)");
539  tmp2D->GetXaxis()->SetTitle("charge (ke-)");
540  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
541  (
542  (Belle2::VxdID)*itSvdSensors,
543  InterHistoAndFill(
544  tmp2D,
545  [this](TH1 * hPtr, const SVDIntercept * inter) {
546 
547  for (auto& it : this->m_SVDClusters)
548  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (it.isUCluster())) {
549  double resid = inter->getCoorU() - it.getPosition(inter->getCoorV());
550  hPtr->Fill(it.getCharge() / 1000., resid);
551  }
552  }
553  )
554  )
555  );
556 
557  name = "hClusterResidV_vs_charge_" + sensorid;
558  title = "Cluster V residual (cm) vs charge " + sensorid;
559  tmp2D = new TH2F(name.c_str(), title.c_str(), 250, 0, 250, 100, -5, 5);
560  tmp2D->GetYaxis()->SetTitle("Cluster V resid (cm)");
561  tmp2D->GetXaxis()->SetTitle("charge (ke-)");
562  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
563  (
564  (Belle2::VxdID)*itSvdSensors,
565  InterHistoAndFill(
566  tmp2D,
567  [this](TH1 * hPtr, const SVDIntercept * inter) {
568 
569  for (auto& it : this->m_SVDClusters)
570  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (!it.isUCluster())) {
571  double resid = inter->getCoorV() - it.getPosition();
572  hPtr->Fill(it.getCharge() / 1000., resid);
573  }
574  }
575  )
576  )
577  );
578 
579  // residual vs time for clusters
580  name = "hClusterResidU_vs_time_" + sensorid;
581  title = "Cluster U residual (cm) vs time " + sensorid;
582  tmp2D = new TH2F(name.c_str(), title.c_str(), 400, -200, 200, 100, -5, 5);
583  tmp2D->GetYaxis()->SetTitle("U resid (cm)");
584  tmp2D->GetXaxis()->SetTitle("time (ns)");
585  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
586  (
587  (Belle2::VxdID)*itSvdSensors,
588  InterHistoAndFill(
589  tmp2D,
590  [this](TH1 * hPtr, const SVDIntercept * inter) {
591 
592  for (auto& it : this->m_SVDClusters)
593  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (it.isUCluster())) {
594  double resid = inter->getCoorU() - it.getPosition(inter->getCoorV());
595  hPtr->Fill(it.getClsTime(), resid);
596  }
597  }
598  )
599  )
600  );
601 
602  name = "hClusterResidV_vs_time_" + sensorid;
603  title = "Cluster V residual (cm) vs time " + sensorid;
604  tmp2D = new TH2F(name.c_str(), title.c_str(), 400, -200, 200, 100, -5, 5);
605  tmp2D->GetYaxis()->SetTitle("Cluster V resid (cm)");
606  tmp2D->GetXaxis()->SetTitle("time (ns)");
607  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
608  (
609  (Belle2::VxdID)*itSvdSensors,
610  InterHistoAndFill(
611  tmp2D,
612  [this](TH1 * hPtr, const SVDIntercept * inter) {
613 
614  for (auto& it : this->m_SVDClusters)
615  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (!it.isUCluster())) {
616  double resid = inter->getCoorV() - it.getPosition();
617  hPtr->Fill(it.getClsTime(), resid);
618  }
619  }
620  )
621  )
622  );
623 
624  // scatter plot: U,V intercept in cm VS U,V cell position
625  name = "hCoorU_vs_UDigit_" + sensorid;
626  title = "U intercept (cm) vs U Digit (ID) " + sensorid;
627  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -5, 5, 1000, -5, 5);
628  tmp2D->GetXaxis()->SetTitle("intercept U coor (cm)");
629  tmp2D->GetYaxis()->SetTitle("digit U coor (cm)");
630  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
631  (
632  (Belle2::VxdID)*itSvdSensors,
633  InterHistoAndFill(
634  tmp2D,
635  [this](TH1 * hPtr, const SVDIntercept * inter) {
636 
637  for (auto& it : this->m_SVDShaperDigits)
638  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (it.isUStrip())) {
639  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
640  hPtr->Fill(inter->getCoorU(), aSensorInfo.getUCellPosition(it.getCellID()));
641  // hPtr->Fill( inter->getCoorU(), it.getVCellID()*75e-4 );
642  }
643  }
644  )
645  )
646  );
647 
648  name = "hCoorV_vs_VDigit_" + sensorid;
649  title = "V intercept (cm) vs V Digit (ID) " + sensorid;
650  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -5, 5, 1000, -5, 5);
651  tmp2D->GetXaxis()->SetTitle("intercept V coor (cm)");
652  tmp2D->GetYaxis()->SetTitle("digi V coor (cm)");
653  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
654  (
655  (Belle2::VxdID)*itSvdSensors,
656  InterHistoAndFill(
657  tmp2D,
658  [this](TH1 * hPtr, const SVDIntercept * inter) {
659 
660  for (auto& it : this->m_SVDShaperDigits) {
661  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (!it.isUStrip())) {
662  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
663  hPtr->Fill(inter->getCoorV(), aSensorInfo.getVCellPosition(it.getCellID()));
664  // hPtr->Fill( inter->getCoorV(), it.getUCellID()*50e-4 );
665  }
666  }
667  }
668  )
669  )
670  );
671 
672 
673 
674 
675 
676 
677  // ------ HISTOGRAMS WITH A FILL PER ROI -------
678  m_ROIDir->cd();
679 
680  // MIN in U and V
681  name = "hminU_" + sensorid;
682  title = "ROI min in U for sensor " + sensorid;
683  hROIDictionary.insert(std::pair< Belle2::VxdID, ROIHistoAndFill >
684  (
685  (Belle2::VxdID)*itSvdSensors,
686  ROIHistoAndFill(
687  new TH1F(name.c_str(), title.c_str(), nPixelsU, 0, nPixelsU),
688  [](TH1 * hPtr, const ROIid * roi) { hPtr->Fill(roi->getMinUid()); }
689  )
690  )
691  );
692  name = "hminV_" + sensorid;
693  title = "ROI min in V for sensor " + sensorid;
694  hROIDictionary.insert(std::pair< Belle2::VxdID, ROIHistoAndFill >
695  (
696  (Belle2::VxdID)*itSvdSensors,
697  ROIHistoAndFill(
698  new TH1F(name.c_str(), title.c_str(), nPixelsV, 0, nPixelsV),
699  [](TH1 * hPtr, const ROIid * roi) { hPtr->Fill(roi->getMinVid()); }
700  )
701  )
702  );
703  //--------------------------
704  // MAX in U and V
705  name = "hmaxU_" + sensorid;
706  title = "ROI max in U for sensor " + sensorid;
707  hROIDictionary.insert(std::pair< Belle2::VxdID, ROIHistoAndFill >
708  (
709  (Belle2::VxdID)*itSvdSensors,
710  ROIHistoAndFill(
711  new TH1F(name.c_str(), title.c_str(), nPixelsU, 0, nPixelsU),
712  [](TH1 * hPtr, const ROIid * roi) { hPtr->Fill(roi->getMaxUid()); }
713  )
714  )
715  );
716  name = "hmaxV_" + sensorid;
717  title = "ROI max in V for sensor " + sensorid;
718  hROIDictionary.insert(std::pair< Belle2::VxdID, ROIHistoAndFill >
719  (
720  (Belle2::VxdID)*itSvdSensors,
721  ROIHistoAndFill(
722  new TH1F(name.c_str(), title.c_str(), nPixelsV, 0, nPixelsV),
723  [](TH1 * hPtr, const ROIid * roi) { hPtr->Fill(roi->getMaxVid()); }
724  )
725  )
726  );
727  //--------------------------
728 
729  // WIDTH in U and V
730  name = "hwidthU_" + sensorid;
731  title = "ROI width in U for sensor " + sensorid;
732  hROIDictionary.insert(std::pair< Belle2::VxdID, ROIHistoAndFill >
733  (
734  (Belle2::VxdID)*itSvdSensors,
735  ROIHistoAndFill(
736  new TH1F(name.c_str(), title.c_str(), nPixelsU, 0, nPixelsU),
737  [](TH1 * hPtr, const ROIid * roi) { hPtr->Fill(roi->getMaxUid() - roi->getMinUid()); }
738  )
739  )
740  );
741  name = "hwidthV_" + sensorid;
742  title = "ROI width in V for sensor " + sensorid;
743  hROIDictionary.insert(std::pair< Belle2::VxdID, ROIHistoAndFill >
744  (
745  (Belle2::VxdID)*itSvdSensors,
746  ROIHistoAndFill(
747  new TH1F(name.c_str(), title.c_str(), nPixelsV, 0, nPixelsV),
748  [](TH1 * hPtr, const ROIid * roi) { hPtr->Fill(roi->getMaxVid() - roi->getMinVid()); }
749  )
750  )
751  );
752 
753  // ROI center
754  name = "hROIcenter_" + sensorid;
755  title = "ROI center " + sensorid;
756  tmp2D = new TH2F(name.c_str(), title.c_str(), nPixelsU, 0, nPixelsU, nPixelsV, 0, nPixelsV);
757  tmp2D->GetXaxis()->SetTitle(" U (ID)");
758  tmp2D->GetYaxis()->SetTitle(" V (ID)");
759  hROIDictionary.insert(std::pair< Belle2::VxdID, ROIHistoAndFill >
760  (
761  (Belle2::VxdID)*itSvdSensors,
762  ROIHistoAndFill(
763  tmp2D,
764  [](TH1 * hPtr, const ROIid * roi) { hPtr->Fill((roi->getMaxUid() + roi->getMinUid()) / 2, (roi->getMaxVid() + roi->getMinVid()) / 2); }
765  )
766  )
767  );
768 
769  //--------------------------
770 
771  ++itSvdSensors;
772  }
773  ++itSvdLadders;
774  }
775  ++itSvdLayers;
776  }
777 
778 }
779 
780 void SVDROIDQMModule::fillSensorInterHistos(const SVDIntercept* inter)
781 {
782 
783  auto its = hInterDictionary.equal_range(inter->getSensorID());
784 
785  for (auto it = its.first; it != its.second; ++it) {
786  InterHistoAndFill aInterHistoAndFill = it->second;
787  aInterHistoAndFill.second(aInterHistoAndFill.first, inter);
788  }
789 
790 }
791 
793 {
794 
795  auto its = hROIDictionary.equal_range(roi->getSensorID());
796 
797  for (auto it = its.first; it != its.second; ++it) {
798  ROIHistoAndFill aROIHistoAndFill = it->second;
799  aROIHistoAndFill.second(aROIHistoAndFill.first, roi);
800  }
801 
802  auto itsEvt = hROIDictionaryEvt.equal_range(roi->getSensorID());
803  for (auto it = itsEvt.first; it != itsEvt.second; ++it)
804  (it->second).accumulate(roi, (it->second).value);
805 }
806 
808 {
809 
810  hCellU->Scale((double)1 / n_events);
811  hCellV->Scale((double)1 / n_events);
812 
813  for (auto it = hROIDictionaryEvt.begin(); it != hROIDictionaryEvt.end(); ++it)
814  delete &(it->second);
815 }
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
ROIid stores the U and V ids and the sensor id of the Region Of Interest.
Definition: ROIid.h:25
SVDIntercept stores the U,V coordinates and uncertainties of the intersection of a track with an SVD ...
Definition: SVDIntercept.h:22
std::unordered_multimap< Belle2::VxdID, ROIHistoAccumulateAndFill &, std::function< size_t(const Belle2::VxdID &) > > hROIDictionaryEvt
map of histograms to be filled once per event
int m_numModules
number of hardware modules
std::unordered_multimap< Belle2::VxdID, InterHistoAndFill, std::function< size_t(const Belle2::VxdID &)> > hInterDictionary
map of histograms to be filled once per intercept
std::pair< TH1 *, std::function< void(TH1 *, const SVDIntercept *) > > InterHistoAndFill
typedef: histograms to be filled once per intercept + filling function
std::string m_InterceptsName
Name of the SVDIntercept StoreArray.
StoreArray< SVDCluster > m_SVDClusters
svd cluster store array
bool m_plotRecoDigits
Produce plots for SVDRecoDigits when True.
VXD::GeoCache & m_geoCache
the geo cache instance
TH1F * hnROIs
number of ROIs
int n_events
number of events
void initialize() override
register histograms
TDirectory * m_ROIDir
ROI directory in the root file.
void fillSensorInterHistos(const SVDIntercept *inter)
fill histograms per sensor, filled once per intercept
void createHistosDictionaries()
create the dictionary
std::string m_SVDRecoDigitsName
reco digit list name
std::string m_SVDClustersName
cluster list name
void fillSensorROIHistos(const ROIid *roi)
fill histograms per sensor, filled once per ROI
void event() override
fill per-event histograms
void endRun() override
fill per-run histograms
SVDROIDQMModule()
Constructor defining the parameters.
int m_specificLayer
specific layer selected for which to produce the plots.
std::pair< TH1 *, std::function< void(TH1 *, const ROIid *) > > ROIHistoAndFill
typedef: histograms to be filled once per roi + filling function
TH1F * hredFactor
reduction factor
StoreArray< SVDShaperDigit > m_SVDShaperDigits
shaper digit store array
StoreArray< SVDRecoDigit > m_SVDRecoDigits
reco digit store array
std::string m_ROIsName
Name of the ROIid StoreArray.
TH1F * hnInter
number of intercpets
TDirectory * m_InterDir
intercepts directory in the root file
std::unordered_multimap< Belle2::VxdID, ROIHistoAndFill, std::function< size_t(const Belle2::VxdID &)> > hROIDictionary
map of histograms to be filled once per roi
TH1F * harea
ROis area.
StoreArray< SVDIntercept > m_Intercepts
SVDIntercept Store Arrays.
std::string m_SVDShaperDigitsName
shaper digit list name
void defineHisto() override
define histograms
StoreArray< ROIid > m_ROIs
ROis store array.
bool isOptional(const std::string &name="")
Tell the DataStore about an optional input.
double getCoorV() const
return the V coordinate of the intercept
Definition: VXDIntercept.h:60
VxdID::baseType getSensorID() const
return the sensor ID
Definition: VXDIntercept.h:66
double getCoorU() const
return the U coordinate of the intercept
Definition: VXDIntercept.h:59
const std::set< Belle2::VxdID > getLayers(SensorInfoBase::SensorType sensortype=SensorInfoBase::VXD)
Return a set of all known Layers.
Definition: GeoCache.cc:176
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:67
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:204
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:193
Base class to provide Sensor Information for PXD and SVD.
int getVCells() const
Return number of pixel/strips in v direction.
int getUCells() const
Return number of pixel/strips in u direction.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
REG_MODULE(arichBtest)
Register the Module.
Abstract base class for different kinds of events.
struct: histograms to be filled once per event + filling fucntion + accumulate function
std::function< void(TH1 *, double &) > fill
fill function