Belle II Software development
SVDDQMInjectionModule.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/svdDQM/SVDDQMInjectionModule.h>
10#include "TDirectory.h"
11
12using namespace std;
13using namespace Belle2;
14using namespace Belle2::SVD;
15using namespace Belle2::VXD;
16
17//-----------------------------------------------------------------
18// Register the Module
19//-----------------------------------------------------------------
20REG_MODULE(SVDDQMInjection);
21
22//-----------------------------------------------------------------
23// Implementation
24//-----------------------------------------------------------------
25
27{
28 //Set module properties
29 setDescription("Monitor SVD Occupancy after Injection.");
31 addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms will be placed.",
32 std::string("SVDInjection"));
33 addParam("ShaperDigits", m_SVDShaperDigitsName, "Name of SVD ShaperDigits to count occupancy - usually ZS5 strips.",
34 std::string(""));
35}
36
38{
39 TDirectory* oldDir = gDirectory;
40 if (m_histogramDirectoryName != "") {
41 oldDir->mkdir(m_histogramDirectoryName.c_str());// do not rely on return value, might be ZERO
42 oldDir->cd(m_histogramDirectoryName.c_str());//changing to the right directory
43 }
44
45 m_hOccAfterInjLER = new TH1F("SVDOccInjLER", "SVDOccInjLER/Time;Time in #mus;Count/Time (5 #mus bins)", 4000, 0, 20000);
46 m_hOccAfterInjHER = new TH1F("SVDOccInjHER", "SVDOccInjHER/Time;Time in #mus;Count/Time (5 #mus bins)", 4000, 0, 20000);
47 m_hTrgOccAfterInjLER = new TH1F("SVDTrgOccInjLER", "SVDTrgOccInjLER/Time;Time in #mus;Triggers/Time (5 #mus bins)", 4000, 0,
48 20000);
49 m_hTrgOccAfterInjHER = new TH1F("SVDTrgOccInjHER", "SVDTrgOccInjHER/Time;Time in #mus;Triggers/Time (5 #mus bins)", 4000, 0,
50 20000);
51 m_hMaxOccAfterInjLER = new TH1F("SVDMaxOccInjLER", "SVDMaxOccInjLER/Time;Time in #mus;Triggers/Time (5 #mus bins)", 4000, 0,
52 20000);
53 m_hMaxOccAfterInjHER = new TH1F("SVDMaxOccInjHER", "SVDMaxOccInjHER/Time;Time in #mus;Triggers/Time (5 #mus bins)", 4000, 0,
54 20000);
55 m_hBunchNumVSNStrips = new TH2F("SVDBunchNumVSNStrips", "SVDBunchNumVSNStrips;Bunch No.;Number of fired strips", 1280, 0, 1280, 10,
56 0,
57 10000);
58
59
60 // cd back to root directory
61 oldDir->cd();
62}
63
65{
66 REG_HISTOGRAM
67
68 m_rawTTD.isOptional();
70}
71
73{
74 // Assume that everything is non-zero ;-)
75 m_hOccAfterInjLER->Reset();
76 m_hOccAfterInjHER->Reset();
77 m_hTrgOccAfterInjLER->Reset();
78 m_hTrgOccAfterInjHER->Reset();
79 m_hMaxOccAfterInjLER->Reset();
80 m_hMaxOccAfterInjHER->Reset();
81 m_hBunchNumVSNStrips->Reset();
82}
83
85{
86
87 if (!m_rawTTD.isValid()) {
88 B2WARNING("Missing RawFTSW, SVDDQMInjection is skipped.");
89 return;
90 }
91
92 if (!m_digits.isValid()) {
93 B2WARNING("Missing " << m_SVDShaperDigitsName << ", SVDDQMInjection is skipped.");
94 return;
95 }
96
97
98 for (auto& it : m_rawTTD) {
99 B2DEBUG(29, "TTD FTSW : " << hex << it.GetTTUtime(0) << " " << it.GetTTCtime(0) << " EvtNr " << it.GetEveNo(0) << " Type " <<
100 (it.GetTTCtimeTRGType(0) & 0xF) << " TimeSincePrev " << it.GetTimeSincePrevTrigger(0) << " TimeSinceInj " <<
101 it.GetTimeSinceLastInjection(0) << " IsHER " << it.GetIsHER(0) << " Bunch " << it.GetBunchNumber(0));
102
103 // get last injection time
104 auto difference = it.GetTimeSinceLastInjection(0);
105 // check time overflow, too long ago
106 if (difference != 0x7FFFFFFF) {
107 unsigned int allV = 0;
108 unsigned int allU = 0;
109 unsigned int nHitsL3v = 0;
110 unsigned int nHitsL3u = 0;
111 for (auto& p : m_digits) {
112 if (p.isUStrip()) allU++;
113 else allV++;
114
115 if (p.getSensorID().getLayerNumber() != 3) continue;
116 if (p.isUStrip())
117 nHitsL3u++;
118 else nHitsL3v++;
119 }
120
121 B2DEBUG(29, "L3 V = " << nHitsL3v << ", L3 U = " << nHitsL3u << ", all V = " << allV << ", all U = " << allU);
122
123 //choose your counter
124 unsigned int counter = nHitsL3u;
125
126 float diff2 = difference / 127.; // 127MHz clock ticks to us, inexact rounding
127 if (it.GetIsHER(0)) {
128 m_hOccAfterInjHER->Fill(diff2, counter);
129 m_hTrgOccAfterInjHER->Fill(diff2);
130 auto bin = m_hMaxOccAfterInjHER->FindBin(diff2);
131 auto value = m_hMaxOccAfterInjHER->GetBinContent(bin);
132 if (counter > value) m_hMaxOccAfterInjHER->SetBinContent(bin, counter);
133 } else {
134 m_hOccAfterInjLER->Fill(diff2, counter);
135 m_hTrgOccAfterInjLER->Fill(diff2);
136 auto bin = m_hMaxOccAfterInjLER->FindBin(diff2);
137 auto value = m_hMaxOccAfterInjLER->GetBinContent(bin);
138 if (counter > value) m_hMaxOccAfterInjLER->SetBinContent(bin, counter);
139
140 }
141 m_hBunchNumVSNStrips->Fill(it.GetBunchNumber(0), allU + allV);
142 }
143
144 break;
145 }
146}
HistoModule()
Constructor.
Definition HistoModule.h:32
void setDescription(const std::string &description)
Sets the description of the module.
Definition Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition Module.h:80
void initialize() override final
initialize function
TH1F * m_hOccAfterInjLER
Histogram Occupancy after LER injection.
TH1F * m_hMaxOccAfterInjLER
Histogram Max Occupancy after LER injection.
TH1F * m_hMaxOccAfterInjHER
Histogram Max Occupancy after HER injection.
StoreArray< RawFTSW > m_rawTTD
Input array for DAQ Status.
void defineHisto() override final
defineHisto function
void event() override final
event function
TH1F * m_hTrgOccAfterInjHER
Histogram for Nr Entries (=Triggrs) for normalization after HER injection.
TH1F * m_hOccAfterInjHER
Histogram Occupancy after HER injection.
std::string m_histogramDirectoryName
Name of the histogram directory in ROOT file.
StoreArray< SVDShaperDigit > m_digits
Input array for SVD Raw Hits.
void beginRun() override final
beginRun function
VXD::GeoCache & m_vxdGeometry
the VXD geometry
TH2F * m_hBunchNumVSNStrips
Histogram Bunch Number VS Number of strips.
SVDDQMInjectionModule()
Constructor defining the parameters.
std::string m_SVDShaperDigitsName
The name of the StoreArray of SVDShaperDigit to be generated.
TH1F * m_hTrgOccAfterInjLER
Histogram for Nr Entries (=Triggrs) for normalization after LER injection.
Class to facilitate easy access to sensor information of the VXD like coordinate transformations or p...
Definition GeoCache.h:39
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition Module.h:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition Module.h:649
Namespace to encapsulate code needed for simulation and reconstrucion of the SVD.
Namespace to provide code needed by both Vertex Detectors, PXD and SVD, and also testbeam telescopes.
Definition GeoCache.h:34
Abstract base class for different kinds of events.
STL namespace.