Belle II Software development
SVDDQMDoseModule.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/SVDDQMDoseModule.h>
10// Includes for the forward declarations
11#include <rawdata/dataobjects/RawFTSW.h>
12#include <svd/dataobjects/SVDShaperDigit.h>
13#include <mdst/dataobjects/TRGSummary.h>
14#include <vxd/dataobjects/VxdID.h>
15// #include <svd/calibration/SVDFADCMaskedStrips.h>
16// Other includes
17#include <vxd/geometry/GeoCache.h>
18#include <svd/geometry/SensorInfo.h>
19#include <TDirectory.h>
20#include <TString.h>
21#include <TMath.h>
22#include <cmath>
23
24using namespace std;
25using namespace Belle2;
26using namespace Belle2::SVD;
27
28REG_MODULE(SVDDQMDose);
29
31{
32 setDescription("The SVD dose-monitoring DQM module. Fills histograms of the SVD's instantaneous occupancy and "
33 "of SVD occupancy vs time since last injection and time in beam revolution cycle.");
35 addParam("eventTypeFilter", m_eventFilter,
36 "Types of events to include in the plots (1 = less than noInjectionTimeout after HER injection, "
37 "2 = less than noInjectionTimeout after LER injection, 4 = more than noInjectionTimeout after any "
38 "injection; bitwise or combinations are possible; see SVDDQMDoseModule::EEventType).", 7U);
39 addParam("histogramDirectoryName", m_histogramDirectoryName,
40 "Name of the directory where histograms will be placed in the ROOT file.",
41 std::string("SVDDose"));
42 addParam("offlineZSShaperDigits", m_SVDShaperDigitsName,
43 "Name of the SVDShaperDigits to use for computing occupancy (default is SVDShaperDigitsZS5).",
44 std::string("SVDShaperDigitsZS5"));
45 addParam("noInjectionTimeout", m_noInjectionTime,
46 "Time (microseconds) since last injection after which an event is considered \"No Injection\". "
47 "Also the limit for the x axis of the 2D histograms.",
48 30e3);
50 addParam("trgTypes", m_trgTypes,
51 "Trigger types for event selection. Empty to select everything. "
52 "Default is only Poisson w/o inj. veto.",
54}
55
57{
58 TDirectory* oldDir = gDirectory;
59 if (m_histogramDirectoryName != "") {
60 oldDir->mkdir(m_histogramDirectoryName.c_str());
61 oldDir->cd(m_histogramDirectoryName.c_str());
62 }
63
64 // Round noInjectionTimeout to the nearest multiple of the revolution cycle
66
67 h_nEvtsVsTime = new TH2F(
68 "SVDEvtsVsTime",
69 "SVD Events;Time since last injection [#mus];Time in beam cycle [#mus];Events / bin",
70 500, 0, m_noInjectionTime, 100, 0, c_revolutionTime);
71
72 m_groupOccupanciesU.reserve(c_sensorGroups.size()); // Allocate memory only once
73 TString name = "SVDInstOccu_";
74 TString title = "SVD Instantaneous Occupancy ";
75 TString axisTitle = ";Occupancy [%];Count / bin";
76 for (const SensorGroup& group : c_sensorGroups) {
77 m_groupOccupanciesU.push_back(
78 new TH1F(name + group.nameSuffix + "U",
79 title + group.titleSuffix + " U-side" + axisTitle,
80 group.nBins, group.xMin, group.xMax));
81 }
82
83 m_groupNHitsU.reserve(c_sensorGroups.size()); // Allocate memory only once
84 name = "SVDHitsVsTime_";
85 title = "SVD Hits ";
86 axisTitle = ";Time since last injection [#mus];Time in beam cycle[#mus];Hits / bin";
87 for (const SensorGroup& group : c_sensorGroups) {
88 m_groupNHitsU.push_back(
89 new TH2F(name + group.nameSuffix + "U",
90 title + group.titleSuffix + " U-side" + axisTitle,
91 500, 0, m_noInjectionTime, 100, 0, c_revolutionTime));
92 }
93
94 // Nbins for 1D histos such that bin width = beam_revolution_time / 2
95 int nb1 = TMath::Nint(m_noInjectionTime * 2.0 / c_revolutionTime);
96 h_nEvtsVsTime1 = new TH1F(
97 "SVDEvtsVsTime1", "SVD Events;Time since last injection [#mus];Events / bin",
98 nb1, 0, m_noInjectionTime);
99
100 m_groupNHits1U.reserve(c_sensorGroups.size()); // Allocate memory only once
101 name = "SVDHitsVsTime1_";
102 title = "SVD Hits ";
103 axisTitle = ";Time since last injection [#mus];Hits / bin";
104 for (const SensorGroup& group : c_sensorGroups) {
105 m_groupNHits1U.push_back(
106 new TH1F(name + group.nameSuffix + "U",
107 title + group.titleSuffix + " U-side" + axisTitle,
108 nb1, 0, m_noInjectionTime));
109 }
110
111 // Include directory name in title: this histogram's canvas is made automatically,
112 // so no analysis modules changes its title to show the event selection used.
113 title = "SVDBunchNumVSNStrips - ";
115 title += ";Bunch No.;Number of fired strips;Events / bin";
116 h_bunchNumVsNHits = new TH2F("SVDBunchNumVSNStrips", title, 1280, 0, 1280, 10, 0, 10000);
117
118 oldDir->cd();
119}
120
122{
123 REG_HISTOGRAM
124
125 // Parameters
126 m_rawTTD.isOptional();
128 m_trgSummary.isOptional();
129
130 // Total number of strips per group
131 static bool nStripsComputed = false; // Compute only once
132 if (nStripsComputed)
133 return;
134 nStripsComputed = true;
136 for (const auto& layer : geo.getLayers(VXD::SensorInfoBase::SVD)) {
137 for (const auto& ladder : geo.getLadders(layer)) {
138 for (const auto& sensor : geo.getSensors(ladder)) {
139 const auto& sInfo = geo.getSensorInfo(sensor);
140 for (const SensorGroup& group : c_sensorGroups) {
141 if (group.contains(sensor)) {
142 // TODO exclude strips that are masked on FADC? It shouldn't matter much...
143 group.nStripsU += sInfo.getUCells();
144 }
145 }
146 }
147 }
148 }
149}
150
152{
153 h_nEvtsVsTime->Reset();
154 h_nEvtsVsTime1->Reset();
155 for (const auto& histPtr : m_groupOccupanciesU)
156 histPtr->Reset();
157 for (const auto& histPtr : m_groupNHitsU)
158 histPtr->Reset();
159 for (const auto& histPtr : m_groupNHits1U)
160 histPtr->Reset();
161 h_bunchNumVsNHits->Reset();
162}
163
165{
166 // Allocate only once, especially good for the vectors
167 // static VXD::GeoCache& geo = VXD::GeoCache::getInstance();
168 static vector<int> groupHitsU(c_sensorGroups.size(), 0);
169
170 if (m_trgTypes.size()) { // If not trg types are given, take everything
171 if (!m_trgSummary.isValid()) {
172 B2WARNING("Missing TRGSummary, SVDDQMDose is skipped.");
173 return;
174 }
175 auto ttyp = m_trgSummary->getTimType();
176 bool discardEvent = true;
177 for (const auto& ttyp2 : m_trgTypes) {
178 if (ttyp == ttyp2) {
179 discardEvent = false;
180 break;
181 }
182 }
183 if (discardEvent)
184 return;
185 }
186
187 if (!m_rawTTD.isValid()) {
188 B2WARNING("Missing RawFTSW, SVDDQMDose is skipped.");
189 return;
190 }
191 if (!m_digits.isValid()) {
192 B2WARNING("Missing SVDShaperDigit " << m_SVDShaperDigitsName
193 << ", SVDDQMDose is skipped.");
194 return;
195 }
196
197 if (m_rawTTD.getEntries() == 0)
198 return;
199 RawFTSW* theTTD = m_rawTTD[0];
200 const double timeSinceInj = theTTD->GetTimeSinceLastInjection(0) / c_globalClock;
201 const bool isHER = theTTD->GetIsHER(0);
202 const EEventType eventType = timeSinceInj > m_noInjectionTime ? c_NoInjection : (isHER ? c_HERInjection : c_LERInjection);
203 if (((unsigned int)eventType & m_eventFilter) == 0U)
204 return;
205 const double timeInCycle = timeSinceInj - (int)(timeSinceInj / c_revolutionTime) * c_revolutionTime;
206
207 // Reset counters
208 for (int& count : groupHitsU) count = 0;
209
210 // Count hits
211 for (const SVDShaperDigit& hit : m_digits) {
212 const VxdID& sensorID = hit.getSensorID();
213 if (hit.isUStrip()) {
214 for (unsigned int i = 0; i < c_sensorGroups.size(); i++) {
215 if (c_sensorGroups[i].contains(sensorID)) {
216 groupHitsU[i]++; // For instantaneous occupancy
217 m_groupNHitsU[i]->Fill(timeSinceInj, timeInCycle);
218 m_groupNHits1U[i]->Fill(timeSinceInj);
219 }
220 }
221 }
222 }
223
224 // Count events
225 h_nEvtsVsTime->Fill(timeSinceInj, timeInCycle);
226 h_nEvtsVsTime1->Fill(timeSinceInj);
227
228 // Compute instantaneous occupancy
229 for (unsigned int i = 0; i < c_sensorGroups.size(); i++)
230 m_groupOccupanciesU[i]->Fill(groupHitsU[i] * 100.0 / c_sensorGroups[i].nStripsU);
231
232 // Bunch num vs fired strips
233 h_bunchNumVsNHits->Fill(theTTD->GetBunchNumber(0), m_digits.getEntries());
234}
235
236const std::vector<SVDDQMDoseModule::SensorGroup> SVDDQMDoseModule::c_sensorGroups = {
237 {"L3XX", "L3", c_defaultNBins, c_defaultOccuMin, c_defaultOccuMax, [](const VxdID & s) { return s.getLayerNumber() == 3; }},
238 {"L4XX", "L4", c_defaultNBins, c_defaultOccuMin, c_defaultOccuMax, [](const VxdID & s) { return s.getLayerNumber() == 4; }},
239 {"L5XX", "L5", c_defaultNBins, c_defaultOccuMin, c_defaultOccuMax, [](const VxdID & s) { return s.getLayerNumber() == 5; }},
240 {"L6XX", "L6", c_defaultNBins, c_defaultOccuMin, c_defaultOccuMax, [](const VxdID & s) { return s.getLayerNumber() == 6; }},
241 {"L31X", "L3.1", c_defaultNBins, c_defaultOccuMin, c_defaultOccuMax, [](const VxdID & s) { return s.getLayerNumber() == 3 && s.getLadderNumber() == 1; }},
242 {"L32X", "L3.2", c_defaultNBins, c_defaultOccuMin, c_defaultOccuMax, [](const VxdID & s) { return s.getLayerNumber() == 3 && s.getLadderNumber() == 2; }}
243};
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
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
The Raw FTSW class.
Definition: RawFTSW.h:30
int GetIsHER(int n)
HER injection = 1 or LER injection = 0.
Definition: RawFTSW.h:173
unsigned int GetTimeSinceLastInjection(int n)
Get time since the last injection.
Definition: RawFTSW.h:180
unsigned int GetBunchNumber(int n)
Get a bunch number.
Definition: RawFTSW.h:194
The SVD ShaperDigit class.
TH1F * h_nEvtsVsTime1
Hist of the total evts in each time bin (1D, time since inj.
unsigned int m_eventFilter
Bitmask for event type selection, see EEventType .
void initialize() override final
Overrides HistoModule::initialize.
TH2F * h_bunchNumVsNHits
Hist of bunch number vs number of fired strips (copied from SVDDQMInjection).
static const std::vector< SensorGroup > c_sensorGroups
List of interesting groups of sensors to average over.
StoreArray< RawFTSW > m_rawTTD
Input: DAQ status.
std::vector< TH1F * > m_groupNHits1U
Hists of the number of hits in each time bin (1D) per sensor group, U-side.
void defineHisto() override final
Overrides HistoModule::defineHisto.
StoreObjPtr< TRGSummary > m_trgSummary
Input: trigger type.
SVDDQMDoseModule()
Default constructor, defines parameters.
static constexpr double c_revolutionTime
Beam revolution time in microseconds (approximated).
void event() override final
Overrides HistoModule::event.
static constexpr double c_globalClock
Approximated global clock frequency in MHz.
double m_noInjectionTime
After this time (in microseconds) from last injection the event falls in the "No Injection" category.
std::vector< int > m_trgTypes
Trigger types to accept (all if the vector is empty).
std::string m_histogramDirectoryName
Name of the histograms' directory in the ROOT file.
std::vector< TH2F * > m_groupNHitsU
Hists of the number of hits in each time bin per sensor group, U-side.
StoreArray< SVDShaperDigit > m_digits
Input: raw hits.
std::vector< TH1F * > m_groupOccupanciesU
Hists of the instantaneous occupancy per sensor group (see c_sensorGroups), U-side.
TH2F * h_nEvtsVsTime
Hist of the total evts in each time bin (time since inj.
EEventType
Bits definition for the bitmask that selects the events to put in the histograms.
void beginRun() override final
Overrides HistoModule::beginRun.
std::string m_SVDShaperDigitsName
Name of the StoreArray of SVDShaperDigit to use.
@ TTYP_POIS
poisson random trigger
Definition: TRGSummary.h:73
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: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
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
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
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
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:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Namespace to encapsulate code needed for simulation and reconstrucion of the SVD.
Definition: GeoSVDCreator.h:23
Abstract base class for different kinds of events.
STL namespace.
A struct to define non-trivial histograms in a human-readable way.