Belle II Software  release-05-01-25
SVDDQMInjectionModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2019 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Giulia Casarosa *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <svd/modules/svdDQM/SVDDQMInjectionModule.h>
12 #include "TDirectory.h"
13 
14 using namespace std;
15 using namespace Belle2;
16 using namespace Belle2::SVD;
17 using namespace Belle2::VXD;
18 
19 //-----------------------------------------------------------------
20 // Register the Module
21 //-----------------------------------------------------------------
22 REG_MODULE(SVDDQMInjection)
23 
24 //-----------------------------------------------------------------
25 // Implementation
26 //-----------------------------------------------------------------
27 
28 SVDDQMInjectionModule::SVDDQMInjectionModule() : HistoModule() , m_vxdGeometry(VXD::GeoCache::getInstance())
29 {
30  //Set module properties
31  setDescription("Monitor SVD Occupancy after Injection");
32  setPropertyFlags(c_ParallelProcessingCertified);
33  addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms will be placed",
34  std::string("SVDInjection"));
35  addParam("ShaperDigits", m_SVDShaperDigitsName, "Name of SVD ShaperDigits to count occupancy", std::string(""));
36 }
37 
38 void SVDDQMInjectionModule::defineHisto()
39 {
40  TDirectory* oldDir = gDirectory;
41  if (m_histogramDirectoryName != "") {
42  oldDir->mkdir(m_histogramDirectoryName.c_str());// do not rely on return value, might be ZERO
43  oldDir->cd(m_histogramDirectoryName.c_str());//changing to the right directory
44  }
45 
46  m_hOccAfterInjLER = new TH1F("SVDOccInjLER", "SVDOccInjLER/Time;Time in #mus;Count/Time (5 #mus bins)", 4000, 0, 20000);
47  m_hOccAfterInjHER = new TH1F("SVDOccInjHER", "SVDOccInjHER/Time;Time in #mus;Count/Time (5 #mus bins)", 4000, 0, 20000);
48  m_hTrgOccAfterInjLER = new TH1F("SVDTrgOccInjLER", "SVDTrgOccInjLER/Time;Time in #mus;Triggers/Time (5 #mus bins)", 4000, 0,
49  20000);
50  m_hTrgOccAfterInjHER = new TH1F("SVDTrgOccInjHER", "SVDTrgOccInjHER/Time;Time in #mus;Triggers/Time (5 #mus bins)", 4000, 0,
51  20000);
52  m_hMaxOccAfterInjLER = new TH1F("SVDMaxOccInjLER", "SVDMaxOccInjLER/Time;Time in #mus;Triggers/Time (5 #mus bins)", 4000, 0,
53  20000);
54  m_hMaxOccAfterInjHER = new TH1F("SVDMaxOccInjHER", "SVDMaxOccInjHER/Time;Time in #mus;Triggers/Time (5 #mus bins)", 4000, 0,
55  20000);
56  m_hBunchNumVSNStrips = new TH2F("SVDBunchNumVSNStrips", "SVDBunchNumVSNStrips;Bunch No.;Number of fired strips", 1280, 0, 1280, 10,
57  0,
58  10000);
59 
60 
61  // cd back to root directory
62  oldDir->cd();
63 }
64 
65 void SVDDQMInjectionModule::initialize()
66 {
67  REG_HISTOGRAM
68 
69  m_rawTTD.isOptional();
70  m_digits.isOptional(m_SVDShaperDigitsName);
71 }
72 
73 void SVDDQMInjectionModule::beginRun()
74 {
75  // Assume that everthing is non-yero ;-)
76  m_hOccAfterInjLER->Reset();
77  m_hOccAfterInjHER->Reset();
78  m_hTrgOccAfterInjLER->Reset();
79  m_hTrgOccAfterInjHER->Reset();
80  m_hMaxOccAfterInjLER->Reset();
81  m_hMaxOccAfterInjHER->Reset();
82  m_hBunchNumVSNStrips->Reset();
83 }
84 
85 void SVDDQMInjectionModule::event()
86 {
87 
88  if (!m_rawTTD.isValid()) {
89  B2WARNING("Missing RawFTSW, SVDDQMInjection is skipped.");
90  return;
91  }
92 
93  if (!m_digits.isValid()) {
94  B2WARNING("Missing " << m_SVDShaperDigitsName << ", SVDDQMInjection is skipped.");
95  return;
96  }
97 
98 
99  for (auto& it : m_rawTTD) {
100  B2DEBUG(29, "TTD FTSW : " << hex << it.GetTTUtime(0) << " " << it.GetTTCtime(0) << " EvtNr " << it.GetEveNo(0) << " Type " <<
101  (it.GetTTCtimeTRGType(0) & 0xF) << " TimeSincePrev " << it.GetTimeSincePrevTrigger(0) << " TimeSinceInj " <<
102  it.GetTimeSinceLastInjection(0) << " IsHER " << it.GetIsHER(0) << " Bunch " << it.GetBunchNumber(0));
103 
104  // get last injection time
105  auto difference = it.GetTimeSinceLastInjection(0);
106  // check time overflow, too long ago
107  if (difference != 0x7FFFFFFF) {
108  unsigned int allV = 0;
109  unsigned int allU = 0;
110  unsigned int nHitsL3v = 0;
111  unsigned int nHitsL3u = 0;
112  for (auto& p : m_digits) {
113  if (p.isUStrip()) allU++;
114  else allV++;
115 
116  if (p.getSensorID().getLayerNumber() != 3) continue;
117  if (p.isUStrip())
118  nHitsL3u++;
119  else nHitsL3v++;
120  }
121 
122  B2DEBUG(29, "L3 V = " << nHitsL3v << ", L3 U = " << nHitsL3u << ", all V = " << allV << ", all U = " << allU);
123 
124  //choose your counter
125  unsigned int counter = nHitsL3u;
126 
127  float diff2 = difference / 127.; // 127MHz clock ticks to us, inexact rounding
128  if (it.GetIsHER(0)) {
129  m_hOccAfterInjHER->Fill(diff2, counter);
130  m_hTrgOccAfterInjHER->Fill(diff2);
131  auto bin = m_hMaxOccAfterInjHER->FindBin(diff2);
132  auto value = m_hMaxOccAfterInjHER->GetBinContent(bin);
133  if (counter > value) m_hMaxOccAfterInjHER->SetBinContent(bin, counter);
134  } else {
135  m_hOccAfterInjLER->Fill(diff2, counter);
136  m_hTrgOccAfterInjLER->Fill(diff2);
137  auto bin = m_hMaxOccAfterInjLER->FindBin(diff2);
138  auto value = m_hMaxOccAfterInjLER->GetBinContent(bin);
139  if (counter > value) m_hMaxOccAfterInjLER->SetBinContent(bin, counter);
140 
141  }
142  m_hBunchNumVSNStrips->Fill(it.GetBunchNumber(0), allU + allV);
143  }
144 
145  break;
146  }
147 }
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::VXD
Namespace to provide code needed by both Vertex Detectors, PXD and SVD, and also testbeam telescopes.
Definition: GeoCache.h:36
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::SVD::SVDDQMInjectionModule
The SVD Occupancy after Injection DQM module.
Definition: SVDDQMInjectionModule.h:45
Belle2::SVD
Namespace to encapsulate code needed for simulation and reconstrucion of the SVD.
Definition: GeoSVDCreator.h:35
Belle2::VXD::GeoCache
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.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