Belle II Software  release-06-02-00
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  , hInterDictionary(172, [](const Belle2::VxdID & vxdid) {return (size_t)vxdid.getID(); })
29 , m_numModules(0)
30 , hnInter(nullptr)
31 , harea(nullptr)
32 , hredFactor(nullptr)
33 , hCellU(nullptr)
34 , hCellV(nullptr)
35 , n_events(0)
36 {
37  //Set module properties
38  setDescription("Monitor of the ROIs creation on HLT");
39  setPropertyFlags(c_ParallelProcessingCertified);
40 
41  addParam("SVDShaperDigitsName", m_SVDShaperDigitsName,
42  "name of the list of SVDShaperDigits", std::string(""));
43  addParam("SVDRecoDigitsName", m_SVDRecoDigitsName,
44  "name of the list of SVDRecoDigits", std::string(""));
45  addParam("SVDClustersName", m_SVDClustersName,
46  "name of the list of SVDClusters", std::string(""));
47 
48  addParam("InterceptsName", m_InterceptsName,
49  "name of the list of interceptions", std::string(""));
50 
51  addParam("specificLayer", m_specificLayer,
52  "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.",
53  m_specificLayer);
54 
55  addParam("plotRecoDigits", m_plotRecoDigits,
56  "Set true to produce the plots for RecoDigits (false by default)", m_plotRecoDigits);
57 
58 }
59 
61 {
62 
63  // Create a separate histogram directory and cd into it.
64  TDirectory* oldDir = gDirectory;
65  TDirectory* roiDir = oldDir->mkdir("SVDROIs");
66  m_InterDir = roiDir->mkdir("intercept");
67 
68  hCellU = new TH1F("hCellU", "CellID U", 769, -0.5, 768.5);
69  hCellU->GetXaxis()->SetTitle("U cell ID");
70  hCellV = new TH1F("hCellV", "CellID V", 769, -0.5, 768.5);
71  hCellV->GetXaxis()->SetTitle("V cell ID");
72 
73  m_InterDir->cd();
74  hnInter = new TH1F("hnInter", "number of intercepts", 100, 0, 100);
75 
77 
78  oldDir->cd();
79 
80 }
81 
83 {
84  REG_HISTOGRAM
85 
89  m_Intercepts.isRequired(m_InterceptsName);
90 
91  n_events = 0;
92 
93 
94 }
95 
97 {
98 
99  n_events++;
100 
101  for (auto& it : m_SVDShaperDigits)
102  if (it.isUStrip())
103  hCellU->Fill(it.getCellID());
104  else
105  hCellV->Fill(it.getCellID());
106 
107  hnInter->Fill(m_Intercepts.getEntries());
108 
109  for (auto& it : m_Intercepts)
111 
112 }
113 
115 {
116 
117  // VXD::GeoCache& aGeometry = VXD::GeoCache::getInstance();
118 
119  std::string name; //name of the histogram
120  std::string title; //title of the histogram
121  TH2F* tmp2D; //temporary 2D histo used to set axis title
122  TH1F* tmp1D; //temporary 1D histo used to set axis title
123 
124  m_numModules = 0;
125 
126  std::set<Belle2::VxdID> svdLayers = m_geoCache.getLayers(VXD::SensorInfoBase::SVD);
127  std::set<Belle2::VxdID>::iterator itSvdLayers = svdLayers.begin();
128 
129  if (m_specificLayer >= 3 && m_specificLayer <= 6) {
130  B2INFO("Producing plots for layer: " << m_specificLayer);
131  svdLayers.clear();
132  svdLayers.insert(Belle2::VxdID(m_specificLayer, 0, 0));
133  itSvdLayers = svdLayers.begin();
134  } else {
135  B2INFO("No specific SVD layer (3,4,5,6) selected (m_specificLayer = " << m_specificLayer <<
136  "). Producing plots for all SVD layers.");
137  }
138 
139  while (itSvdLayers != svdLayers.end()) {
140 
141  std::set<Belle2::VxdID> svdLadders = m_geoCache.getLadders(*itSvdLayers);
142  std::set<Belle2::VxdID>::iterator itSvdLadders = svdLadders.begin();
143 
144  while (itSvdLadders != svdLadders.end()) {
145 
146  std::set<Belle2::VxdID> svdSensors = m_geoCache.getSensors(*itSvdLadders);
147  std::set<Belle2::VxdID>::iterator itSvdSensors = svdSensors.begin();
148 
149  while (itSvdSensors != svdSensors.end()) {
150 
151  m_numModules++; //counting the total number of modules
152 
153  std::string sensorid = std::to_string(itSvdSensors->getLayerNumber()) + "_" + std::to_string(
154  itSvdSensors->getLadderNumber()) + "_" +
155  std::to_string(itSvdSensors->getSensorNumber());
156 
157 
158  // ------ HISTOGRAMS WITH A FILL PER INTERCEPT -------
159  m_InterDir->cd();
160 
161  // coor U and V
162  name = "hCoorU_" + sensorid;
163  title = "U coordinate of the extrapolation in U for sensor " + sensorid;
164  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
165  (
166  (Belle2::VxdID)*itSvdSensors,
168  new TH1F(name.c_str(), title.c_str(), 100, -5, 5),
169  [](TH1 * hPtr, const SVDIntercept * inter) { hPtr->Fill(inter->getCoorU()); }
170  )
171  )
172  );
173 
174  name = "hCoorV_" + sensorid;
175  title = "V coordinate of the extrapolation in V for sensor " + sensorid;
176  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
177  (
178  (Belle2::VxdID)*itSvdSensors,
180  new TH1F(name.c_str(), title.c_str(), 100, -5, 5),
181  [](TH1 * hPtr, const SVDIntercept * inter) { hPtr->Fill(inter->getCoorV()); }
182  )
183  )
184  );
185 
186  // Intercept U vs V coordinate
187  name = "hCoorU_vs_CoorV_" + sensorid;
188  title = "U vs V intercept (cm) " + sensorid;
189  tmp2D = new TH2F(name.c_str(), title.c_str(), 100, -5, 5, 100, -5, 5);
190  tmp2D->GetXaxis()->SetTitle("intercept U coor (cm)");
191  tmp2D->GetYaxis()->SetTitle("intercept V coor (cm)");
192  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
193  (
194  (Belle2::VxdID)*itSvdSensors,
196  tmp2D,
197  [](TH1 * hPtr, const SVDIntercept * inter) { hPtr->Fill(inter->getCoorU(), inter->getCoorV()); }
198  )
199  )
200  );
201 
202 
203  // sigma U and V
204  name = "hStatErrU_" + sensorid;
205  title = "stat error of the extrapolation in U for sensor " + sensorid;
206  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
207  (
208  (Belle2::VxdID)*itSvdSensors,
210  new TH1F(name.c_str(), title.c_str(), 100, 0, 0.35),
211  [](TH1 * hPtr, const SVDIntercept * inter) { hPtr->Fill(inter->getSigmaU()); }
212  )
213  )
214  );
215  name = "hStatErrV_" + sensorid;
216  title = "stat error of the extrapolation in V for sensor " + sensorid;
217  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
218  (
219  (Belle2::VxdID)*itSvdSensors,
221  new TH1F(name.c_str(), title.c_str(), 100, 0, 0.35),
222  [](TH1 * hPtr, const SVDIntercept * inter) { hPtr->Fill(inter->getSigmaV()); }
223  )
224  )
225  );
226 
227  //1D residuals
228  name = "hResidU_" + sensorid;
229  title = "U residuals = intercept - digit, for sensor " + sensorid;
230  tmp1D = new TH1F(name.c_str(), title.c_str(), 1000, -5, 5);
231  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
232  (
233  (Belle2::VxdID)*itSvdSensors,
235  tmp1D,
236  [this](TH1 * hPtr, const SVDIntercept * inter) {
237 
238  for (auto& it : this->m_SVDShaperDigits)
239  if ((int)it.getSensorID() == (int)inter->getSensorID()) {
240  if (it.isUStrip()) {
241  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
242  hPtr->Fill(inter->getCoorU() - aSensorInfo.getUCellPosition(it.getCellID()));
243  }
244  }
245  }
246  )
247  )
248  );
249 
250  name = "hResidV_" + sensorid;
251  title = "V residuals = intercept - digit, for sensor " + sensorid;
252  tmp1D = new TH1F(name.c_str(), title.c_str(), 1000, -5, 5);
253  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
254  (
255  (Belle2::VxdID)*itSvdSensors,
257  tmp1D,
258  [this](TH1 * hPtr, const SVDIntercept * inter) {
259 
260  for (auto& it : this->m_SVDShaperDigits)
261  if ((int)it.getSensorID() == (int)inter->getSensorID()) {
262  if (!it.isUStrip()) {
263  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
264  hPtr->Fill(inter->getCoorV() - aSensorInfo.getVCellPosition(it.getCellID()));
265  }
266  }
267  }
268  )
269  )
270  );
271 
272 
273  //residual U,V vs coordinate U,V
274  name = "hResidU_vs_CoorU_" + sensorid;
275  title = "U residual (cm) vs coor U (cm) " + sensorid;
276  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -5, 5, 1000, -5, 5);
277  tmp2D->GetYaxis()->SetTitle("U resid (cm)");
278  tmp2D->GetXaxis()->SetTitle("U coor (cm)");
279  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
280  (
281  (Belle2::VxdID)*itSvdSensors,
283  tmp2D,
284  [this](TH1 * hPtr, const SVDIntercept * inter) {
285 
286  for (auto& it : this->m_SVDShaperDigits)
287  if (((int)it.getSensorID() == (int)inter->getSensorID()) && it.isUStrip()) {
288  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
289  double resid = inter->getCoorU() - aSensorInfo.getUCellPosition(it.getCellID());
290  hPtr->Fill(inter->getCoorU(), resid);
291  }
292  }
293  )
294  )
295  );
296 
297  name = "hResidV_vs_CoorV_" + sensorid;
298  title = "V residual (cm) vs coor V (cm) " + sensorid;
299  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -5, 5, 1000, -5, 5);
300  tmp2D->GetYaxis()->SetTitle("V resid (cm)");
301  tmp2D->GetXaxis()->SetTitle("V coor (cm)");
302  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
303  (
304  (Belle2::VxdID)*itSvdSensors,
306  tmp2D,
307  [this](TH1 * hPtr, const SVDIntercept * inter) {
308 
309  for (auto& it : this->m_SVDShaperDigits)
310  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (!it.isUStrip())) {
311  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
312  double resid = inter->getCoorV() - aSensorInfo.getVCellPosition(it.getCellID());
313  hPtr->Fill(inter->getCoorV(), resid);
314  }
315  }
316  )
317  )
318  );
319 
320 
321  // RecoDigits
322  //residual vs charge
323  if (m_plotRecoDigits) {
324  name = "hResidU_vs_charge_" + sensorid;
325  title = "U residual (cm) vs charge " + sensorid;
326  tmp2D = new TH2F(name.c_str(), title.c_str(), 250, 0, 250, 100, -5, 5);
327  tmp2D->GetYaxis()->SetTitle("U resid (cm)");
328  tmp2D->GetXaxis()->SetTitle("charge");
329  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
330  (
331  (Belle2::VxdID)*itSvdSensors,
333  tmp2D,
334  [this](TH1 * hPtr, const SVDIntercept * inter) {
335 
336  for (auto& it : this->m_SVDRecoDigits)
337  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (it.isUStrip())) {
338  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
339  double resid = inter->getCoorU() - aSensorInfo.getUCellPosition(it.getCellID());
340  hPtr->Fill(it.getCharge(), resid);
341  }
342  }
343  )
344  )
345  );
346 
347  name = "hResidV_vs_charge_" + sensorid;
348  title = "V residual (cm) vs charge " + sensorid;
349  tmp2D = new TH2F(name.c_str(), title.c_str(), 250, 0, 250, 100, -5, 5);
350  tmp2D->GetYaxis()->SetTitle("V resid (cm)");
351  tmp2D->GetXaxis()->SetTitle("charge");
352  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
353  (
354  (Belle2::VxdID)*itSvdSensors,
356  tmp2D,
357  [this](TH1 * hPtr, const SVDIntercept * inter) {
358 
359  for (auto& it : this->m_SVDRecoDigits)
360  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (!it.isUStrip())) {
361  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
362  double resid = inter->getCoorV() - aSensorInfo.getVCellPosition(it.getCellID());
363  hPtr->Fill(it.getCharge(), resid);
364  }
365  }
366  )
367  )
368  );
369  }
370 
371  // 1D residual for clusters
372  name = "hClusterResidU_" + sensorid;
373  title = "Cluster U residuals = intercept - cluster, for sensor " + sensorid;
374  tmp1D = new TH1F(name.c_str(), title.c_str(), 1000, -5, 5);
375  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
376  (
377  (Belle2::VxdID)*itSvdSensors,
379  tmp1D,
380  [this](TH1 * hPtr, const SVDIntercept * inter) {
381 
382  for (auto& it : this->m_SVDClusters)
383  if ((int)it.getSensorID() == (int)inter->getSensorID()) {
384  if (it.isUCluster()) {
385  hPtr->Fill(inter->getCoorU() - it.getPosition(inter->getCoorV()));
386  }
387  }
388  }
389  )
390  )
391  );
392 
393  name = "hClusterResidV_" + sensorid;
394  title = "Cluster V residuals = intercept - cluster, for sensor " + sensorid;
395  tmp1D = new TH1F(name.c_str(), title.c_str(), 1000, -5, 5);
396  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
397  (
398  (Belle2::VxdID)*itSvdSensors,
400  tmp1D,
401  [this](TH1 * hPtr, const SVDIntercept * inter) {
402 
403  for (auto& it : this->m_SVDClusters)
404  if ((int)it.getSensorID() == (int)inter->getSensorID()) {
405  if (!it.isUCluster()) {
406  hPtr->Fill(inter->getCoorV() - it.getPosition());
407  }
408  }
409  }
410  )
411  )
412  );
413 
414  //residual U,V vs coordinate U,V for clusters
415  name = "hClusterResidU_vs_CoorU_" + sensorid;
416  title = "Cluster U residual (cm) vs coor U (cm) " + sensorid;
417  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -5, 5, 1000, -5, 5);
418  tmp2D->GetYaxis()->SetTitle("Cluster U resid (cm)");
419  tmp2D->GetXaxis()->SetTitle("U coor (cm)");
420  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
421  (
422  (Belle2::VxdID)*itSvdSensors,
423  InterHistoAndFill(
424  tmp2D,
425  [this](TH1 * hPtr, const SVDIntercept * inter) {
426 
427  for (auto& it : this->m_SVDClusters)
428  if (((int)it.getSensorID() == (int)inter->getSensorID()) && it.isUCluster()) {
429  double resid = inter->getCoorU() - it.getPosition(inter->getCoorV());
430  hPtr->Fill(inter->getCoorU(), resid);
431  }
432  }
433  )
434  )
435  );
436 
437  name = "hClusterResidV_vs_CoorV_" + sensorid;
438  title = "Cluster V residual (cm) vs coor V (cm) " + sensorid;
439  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -5, 5, 1000, -5, 5);
440  tmp2D->GetYaxis()->SetTitle("Cluster V resid (cm)");
441  tmp2D->GetXaxis()->SetTitle("V coor (cm)");
442  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
443  (
444  (Belle2::VxdID)*itSvdSensors,
445  InterHistoAndFill(
446  tmp2D,
447  [this](TH1 * hPtr, const SVDIntercept * inter) {
448 
449  for (auto& it : this->m_SVDClusters)
450  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (!it.isUCluster())) {
451  double resid = inter->getCoorV() - it.getPosition();
452  hPtr->Fill(inter->getCoorV(), resid);
453  }
454  }
455  )
456  )
457  );
458 
459 
460  //residual vs charge for clusters
461  name = "hClusterResidU_vs_charge_" + sensorid;
462  title = "Cluster U residual (cm) vs charge " + sensorid;
463  tmp2D = new TH2F(name.c_str(), title.c_str(), 250, 0, 250, 100, -5, 5);
464  tmp2D->GetYaxis()->SetTitle("U resid (cm)");
465  tmp2D->GetXaxis()->SetTitle("charge (ke-)");
466  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
467  (
468  (Belle2::VxdID)*itSvdSensors,
469  InterHistoAndFill(
470  tmp2D,
471  [this](TH1 * hPtr, const SVDIntercept * inter) {
472 
473  for (auto& it : this->m_SVDClusters)
474  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (it.isUCluster())) {
475  double resid = inter->getCoorU() - it.getPosition(inter->getCoorV());
476  hPtr->Fill(it.getCharge() / 1000., resid);
477  }
478  }
479  )
480  )
481  );
482 
483  name = "hClusterResidV_vs_charge_" + sensorid;
484  title = "Cluster V residual (cm) vs charge " + sensorid;
485  tmp2D = new TH2F(name.c_str(), title.c_str(), 250, 0, 250, 100, -5, 5);
486  tmp2D->GetYaxis()->SetTitle("Cluster V resid (cm)");
487  tmp2D->GetXaxis()->SetTitle("charge (ke-)");
488  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
489  (
490  (Belle2::VxdID)*itSvdSensors,
491  InterHistoAndFill(
492  tmp2D,
493  [this](TH1 * hPtr, const SVDIntercept * inter) {
494 
495  for (auto& it : this->m_SVDClusters)
496  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (!it.isUCluster())) {
497  double resid = inter->getCoorV() - it.getPosition();
498  hPtr->Fill(it.getCharge() / 1000., resid);
499  }
500  }
501  )
502  )
503  );
504 
505 
506  // scatter plot: U,V intercept in cm VS U,V cell position
507  name = "hCoorU_vs_UDigit_" + sensorid;
508  title = "U intercept (cm) vs U Digit (ID) " + sensorid;
509  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -5, 5, 1000, -5, 5);
510  tmp2D->GetXaxis()->SetTitle("intercept U coor (cm)");
511  tmp2D->GetYaxis()->SetTitle("digit U coor (cm)");
512  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
513  (
514  (Belle2::VxdID)*itSvdSensors,
515  InterHistoAndFill(
516  tmp2D,
517  [this](TH1 * hPtr, const SVDIntercept * inter) {
518 
519  for (auto& it : this->m_SVDShaperDigits)
520  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (it.isUStrip())) {
521  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
522  hPtr->Fill(inter->getCoorU(), aSensorInfo.getUCellPosition(it.getCellID()));
523  // hPtr->Fill( inter->getCoorU(), it.getVCellID()*75e-4 );
524  }
525  }
526  )
527  )
528  );
529 
530  name = "hCoorV_vs_VDigit_" + sensorid;
531  title = "V intercept (cm) vs V Digit (ID) " + sensorid;
532  tmp2D = new TH2F(name.c_str(), title.c_str(), 1000, -5, 5, 1000, -5, 5);
533  tmp2D->GetXaxis()->SetTitle("intercept V coor (cm)");
534  tmp2D->GetYaxis()->SetTitle("digi V coor (cm)");
535  hInterDictionary.insert(std::pair< Belle2::VxdID, InterHistoAndFill >
536  (
537  (Belle2::VxdID)*itSvdSensors,
538  InterHistoAndFill(
539  tmp2D,
540  [this](TH1 * hPtr, const SVDIntercept * inter) {
541 
542  for (auto& it : this->m_SVDShaperDigits) {
543  if (((int)it.getSensorID() == (int)inter->getSensorID()) && (!it.isUStrip())) {
544  const VXD::SensorInfoBase& aSensorInfo = m_geoCache.getSensorInfo(it.getSensorID());
545  hPtr->Fill(inter->getCoorV(), aSensorInfo.getVCellPosition(it.getCellID()));
546  // hPtr->Fill( inter->getCoorV(), it.getUCellID()*50e-4 );
547  }
548  }
549  }
550  )
551  )
552  );
553 
554 
555  //--------------------------
556 
557  ++itSvdSensors;
558  }
559  ++itSvdLadders;
560  }
561  ++itSvdLayers;
562  }
563 
564 }
565 
566 void SVDROIDQMModule::fillSensorInterHistos(const SVDIntercept* inter)
567 {
568 
569  auto its = hInterDictionary.equal_range(inter->getSensorID());
570 
571  for (auto it = its.first; it != its.second; ++it) {
572  InterHistoAndFill aInterHistoAndFill = it->second;
573  aInterHistoAndFill.second(aInterHistoAndFill.first, inter);
574  }
575 
576 }
577 
578 
580 {
581 
582  hCellU->Scale((double)1 / n_events);
583  hCellV->Scale((double)1 / n_events);
584 }
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
SVDIntercept stores the U,V coordinates and uncertainties of the intersection of a track with an SVD ...
Definition: SVDIntercept.h:22
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
int n_events
number of events
void initialize() override
register histograms
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 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.
StoreArray< SVDShaperDigit > m_SVDShaperDigits
shaper digit store array
StoreArray< SVDRecoDigit > m_SVDRecoDigits
reco digit store array
TH1F * hnInter
number of intercpets
TDirectory * m_InterDir
intercepts directory in the root file
StoreArray< SVDIntercept > m_Intercepts
SVDIntercept Store Arrays.
std::string m_SVDShaperDigitsName
shaper digit list name
void defineHisto() override
define histograms
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: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
Abstract base class for different kinds of events.