Belle II Software development
SVDOccupancyAnalysisModule.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/svdPerformance/SVDOccupancyAnalysisModule.h>
10#include <hlt/softwaretrigger/core/FinalTriggerDecisionCalculator.h>
11
12#include <TMath.h>
13
14using namespace std;
15using namespace Belle2;
16using namespace SoftwareTrigger;
17
18REG_MODULE(SVDOccupancyAnalysis);
19
21{
22
23 setDescription("This module check performances of SVD reconstruction of VXD TB data");
24
25 addParam("outputFileName", m_rootFileName, "Name of output root file.", std::string("SVDOccupancyAnalysis_output.root"));
26
27 addParam("skipHLTRejectedEvents", m_skipRejectedEvents, "If TRUE skip events rejected by HLT", bool(false));
28 addParam("groupNevents", m_group, "Number of events to group", float(10000));
29 addParam("FADCmode", m_FADCmode,
30 "FADC mode: if true the approximation to integer is done", bool(false));
31 addParam("minZScut", m_minZS, "Minimum ZS cut", float(3));
32 addParam("maxZScut", m_maxZS, "Maximum ZS cut", float(6));
33 addParam("pointsZScut", m_pointsZS, "Number of ZS cuts", int(8));
34
35 addParam("ShaperDigitsName", m_ShaperDigitName, "Name of ShaperDigit Store Array.", std::string(""));
36}
37
39{
40
41 m_eventMetaData.isRequired();
43 B2INFO(" ShaperDigits: " << m_ShaperDigitName);
44
45 m_rootFilePtr = new TFile(m_rootFileName.c_str(), "RECREATE");
46
47 m_nEvents = 0;
48
49}
50
51
53{
54
55
56 m_occ_L3U = new TH1F("occL3U", "Occupancy Distribution for L3 U side", m_distr_Nbins, m_distr_min, m_distr_max);
57 m_occ_L3U->GetXaxis()->SetTitle("occupancy(%)");
58 m_occ_L3V = new TH1F("occL3V", "Occupancy Distribution for L3 V side", m_distr_Nbins, m_distr_min, m_distr_max);
59 m_occ_L3V->GetXaxis()->SetTitle("occupancy(%)");
60 m_occ_L4U = new TH1F("occL4U", "Occupancy Distribution for L4 U side", m_distr_Nbins, m_distr_min, m_distr_max);
61 m_occ_L4U->GetXaxis()->SetTitle("occupancy(%)");
62 m_occ_L4V = new TH1F("occL4V", "Occupancy Distribution for L4 V side", m_distr_Nbins, m_distr_min, m_distr_max);
63 m_occ_L4V->GetXaxis()->SetTitle("occupancy(%)");
64 m_occ_L5U = new TH1F("occL5U", "Occupancy Distribution for L5 U side", m_distr_Nbins, m_distr_min, m_distr_max);
65 m_occ_L5U->GetXaxis()->SetTitle("occupancy(%)");
66 m_occ_L5V = new TH1F("occL5V", "Occupancy Distribution for L5 V side", m_distr_Nbins, m_distr_min, m_distr_max);
67 m_occ_L5V->GetXaxis()->SetTitle("occupancy(%)");
68 m_occ_L6U = new TH1F("occL6U", "Occupancy Distribution for L6 U side", m_distr_Nbins, m_distr_min, m_distr_max);
69 m_occ_L6U->GetXaxis()->SetTitle("occupancy(%)");
70 m_occ_L6V = new TH1F("occL6V", "Occupancy Distribution for L6 V side", m_distr_Nbins, m_distr_min, m_distr_max);
71 m_occ_L6V->GetXaxis()->SetTitle("occupancy(%)");
72
73
75
76 //collect the list of all SVD Modules in the geometry here
77 std::vector<VxdID> sensors = geo.getListOfSensors();
78 for (VxdID& aVxdID : sensors) {
79 VXD::SensorInfoBase info = geo.getSensorInfo(aVxdID);
80 if (info.getType() != VXD::SensorInfoBase::SVD) continue;
81 m_SVDModules.push_back(aVxdID); // reorder, sort would be better
82 }
83 std::sort(m_SVDModules.begin(), m_SVDModules.end()); // back to natural order
84
85 m_hit = new SVDSummaryPlots("hits@view", "Number of hits on @view/@side Side");
86
87 TH1F h_dist("dist_L@layerL@ladderS@sensor@view",
88 "Occupancy Distribution (layer @layer, ladder @ladder, sensor @sensor, side@view/@side)",
90 h_dist.GetXaxis()->SetTitle("occupancy (%)");
92
93 TH1F h_occ_768("occ768_L@layerL@ladderS@sensor@view", "Occupancy (layer @layer, ladder @ladder, sensor @sensor, side@view/@side)",
94 768, 0, 768);
95 h_occ_768.GetXaxis()->SetTitle("cellID");
96 TH1F h_occ_512("occ512_L@layerL@ladderS@sensor@view", "Occupancy (layer @layer, ladder @ladder, sensor @sensor, side@view/@side)",
97 512, 0, 512);
98 h_occ_512.GetXaxis()->SetTitle("cellID");
99 m_histo_occ = new SVDHistograms<TH1F>(h_occ_768, h_occ_768, h_occ_768, h_occ_512);
100
101
102 TH1F h_zsVSocc("occVSzs_L@layerL@ladderS@sensor@view",
103 "Average Occupancy VS Zero Suppression (layer @layer, ladder @ladder, sensor @sensor, side@view/@side)", m_pointsZS, m_minZS,
104 m_maxZS);
105 h_zsVSocc.GetXaxis()->SetTitle("ZS cut");
106 m_histo_zsOcc = new SVDHistograms<TH1F>(h_zsVSocc);
107
108
109 TH1F h_zsVSoccSQ("zsVSoccSQ_L@layerL@ladderS@sensor@view",
110 "Average Occupancy VS (ZS cut)^2 (layer @layer, ladder @ladder, sensor @sensor, side@view/@side)", 100, TMath::Power(m_minZS,
111 2) - 5, TMath::Power(m_maxZS, 2));
112 h_zsVSoccSQ.GetXaxis()->SetTitle("(ZS cut)^2");
113 m_histo_zsOccSQ = new SVDHistograms<TH1F>(h_zsVSoccSQ);
114
115
116 TH2F h_occtdep_768("occ768VSevt_L@layerL@ladderS@sensor@view",
117 "Average Occupancy VS Event Number (layer @layer, ladder @ladder, sensor @sensor, side@view/@side)",
118 1000, 0, 1000, 768, 0, 768);
119 h_occtdep_768.GetXaxis()->SetTitle(Form("evt number/%1.0f", m_group));
120 h_occtdep_768.GetYaxis()->SetTitle("cellID");
121
122 TH2F h_occtdep_512("occ512VSevt_L@layerL@ladderS@sensor@view",
123 "Average Occupancy VS Event Number (layer @layer, ladder @ladder, sensor @sensor, side@view/@side)",
124 1000, 0, 1000, 512, 0, 512);
125 h_occtdep_512.GetXaxis()->SetTitle(Form("evt number/%1.0f", m_group));
126 h_occtdep_512.GetYaxis()->SetTitle("cellID");
127
128 m_histo_occtdep = new SVDHistograms<TH2F>(h_occtdep_768, h_occtdep_768, h_occtdep_768, h_occtdep_512);
129
130}
131
133{
134
137 if (!eventAccepted) return;
138 }
139
140 m_nEvents++;
141 int nEvent = m_eventMetaData->getEvent();
142
143 //shaper digits
144 for (int digi = 0 ; digi < m_svdShapers.getEntries(); digi++) {
145
146
147 VxdID::baseType theVxdID = (VxdID::baseType)m_svdShapers[digi]->getSensorID();
148 int side = m_svdShapers[digi]->isUStrip();
149
150 //fill standard occupancy plot, for default zero suppression
151 m_histo_occtdep->fill(theVxdID, side, nEvent / m_group, m_svdShapers[digi]->getCellID());
152
153 m_hit->fill(theVxdID, side, 1);
154
155 m_histo_occ->fill(theVxdID, side, m_svdShapers[digi]->getCellID());
156
157 float noise = m_NoiseCal.getNoise(theVxdID, side, m_svdShapers[digi]->getCellID());
158 float step = (m_maxZS - m_minZS) / m_pointsZS;
159
160 for (int z = 0; z <= m_pointsZS; z++) {
161 int nOKSamples = 0;
162 float cutMinSignal = (m_minZS + step * z) * noise;
163
164 if (m_FADCmode) {
165 cutMinSignal = cutMinSignal + 0.5;
166 cutMinSignal = (int)cutMinSignal;
167 }
168
169
170 Belle2::SVDShaperDigit::APVFloatSamples samples_vec = m_svdShapers[digi]->getSamples();
171
172 for (int k = 0; k < 6; k ++)
173 if (samples_vec[k] > cutMinSignal)
174 nOKSamples++;
175
176 if (nOKSamples > 0) {
177 m_histo_zsOcc->fill(theVxdID, side, m_minZS + z * step);
178 m_histo_zsOccSQ->fill(theVxdID, side, TMath::Power(m_minZS + z * step, 2));
179 }
180 }
181
182 }
183
184 //loop on sensors, fill and clear
185 for (unsigned int i = 0; i < m_SVDModules.size(); i++) {
186 B2DEBUG(10, "module " << i << "," << m_SVDModules[i]);
187 float nStripsV = 512;
188 if (m_SVDModules[i].getLayerNumber() == 3)
189 nStripsV = 768;
190
191 double occU = 100. * m_hit->getValue(m_SVDModules[i], 1) / 768;
192 double occV = 100. * m_hit->getValue(m_SVDModules[i], 0) / nStripsV;
193
194 m_histo_dist->fill(m_SVDModules[i], 1, occU);
195 m_histo_dist->fill(m_SVDModules[i], 0, occV);
196
197 if (m_SVDModules[i].getLayerNumber() == 3) {
198 m_occ_L3U->Fill(occU);
199 m_occ_L3V->Fill(occV);
200 }
201 if (m_SVDModules[i].getLayerNumber() == 4) {
202 m_occ_L4U->Fill(occU);
203 m_occ_L4V->Fill(occV);
204 }
205 if (m_SVDModules[i].getLayerNumber() == 5) {
206 m_occ_L5U->Fill(occU);
207 m_occ_L5V->Fill(occV);
208 }
209 if (m_SVDModules[i].getLayerNumber() == 6) {
210 m_occ_L6U->Fill(occU);
211 m_occ_L6V->Fill(occV);
212 }
213
214 }
215
216 (m_hit->getHistogram(0))->Reset();
217 (m_hit->getHistogram(1))->Reset();
218
219}
220
221
223{
224
225 if (m_rootFilePtr != nullptr) {
226 m_rootFilePtr->cd();
227
228 TDirectory* oldDir = gDirectory;
229
230 m_occ_L3U->Write();
231 m_occ_L3V->Write();
232 m_occ_L4U->Write();
233 m_occ_L4V->Write();
234 m_occ_L5U->Write();
235 m_occ_L5V->Write();
236 m_occ_L6U->Write();
237 m_occ_L6V->Write();
238
240
241 for (auto layer : geoCache.getLayers(VXD::SensorInfoBase::SVD)) {
242 TString layerName = Form("occupancyL%d", layer.getLayerNumber());
243 TDirectory* dir_layer = oldDir->mkdir(layerName.Data());
244 dir_layer->cd();
245 for (auto ladder : geoCache.getLadders(layer))
246 for (Belle2::VxdID sensor : geoCache.getSensors(ladder))
247 for (int view = SVDHistograms<TH1F>::VIndex ; view < SVDHistograms<TH1F>::UIndex + 1; view++) {
248 (m_histo_dist->getHistogram(sensor, view))->Write();
249
250 (m_histo_occ->getHistogram(sensor, view))->Scale(1. / m_nEvents);
251 (m_histo_occ->getHistogram(sensor, view))->Write();
252
253 int nStrips = 768;
254 if (sensor.getLayerNumber() != 3 && view == SVDHistograms<TH1F>::VIndex)
255 nStrips = 512;
256
257 (m_histo_zsOcc->getHistogram(sensor, view))->Scale(1. / m_nEvents / nStrips);
258 (m_histo_zsOcc->getHistogram(sensor, view))->Write();
259 (m_histo_zsOccSQ->getHistogram(sensor, view))->Scale(1. / m_nEvents / nStrips);
260 (m_histo_zsOccSQ->getHistogram(sensor, view))->Write();
261
262 (m_histo_occtdep->getHistogram(sensor, view))->Scale(1. / m_group);
263 (m_histo_occtdep->getHistogram(sensor, view))->Write();
264 }
265 }
266
267 m_rootFilePtr->Close();
268
269 }
270}
271
272
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
template class for SVd histograms
Definition: SVDHistograms.h:24
void fill(const VxdID &vxdID, int view, Types ... args)
fill the histogram for
Definition: SVDHistograms.h:77
H * getHistogram(const VxdID &vxdID, int view)
get a reference to the histogram for
Definition: SVDHistograms.h:56
float getNoise(const VxdID &sensorID, const bool &isU, const unsigned short &strip) const
This is the method for getting the noise.
SVDSummaryPlots * m_hit
hit number summary histogram
int m_pointsZS
number of steps for different ZS cuts
bool m_skipRejectedEvents
if true skip events rejected by HLT
bool m_FADCmode
if true, ZS done with same algorithm as on FADC
SVDHistograms< TH1F > * m_histo_occ
occupancy histograms
std::string m_ShaperDigitName
ShaperDigit StoreArray name.
virtual void initialize() override
check StoreArrays & create rootfile
TH1F * m_occ_L4V
occupancy distribution for L4 V-side sensors
virtual void event() override
fill histograms
SVDHistograms< TH2F > * m_histo_occtdep
occupancy VS event number
SVDHistograms< TH1F > * m_histo_zsOccSQ
occupancy VS ZS cut swuared histograms
std::vector< VxdID > m_SVDModules
IDs of all SVD Modules to iterate over.
SVDNoiseCalibrations m_NoiseCal
SVDNoise calibrations db object.
virtual void endRun() override
write histogrmas
float m_group
number of events to compute occupancy for occ VS time
TH1F * m_occ_L6U
occupancy distribution for L6 U-side sensors
float m_minZS
minimum zero suppresion cut
TH1F * m_occ_L5U
occupancy distribution for L5 U-side sensors
StoreObjPtr< EventMetaData > m_eventMetaData
Event Meta Data StoreObjectPointer.
SVDHistograms< TH1F > * m_histo_dist
occupancy distribution histograms
float m_distr_Nbins
number of bins of occupancy distributions
virtual void beginRun() override
create histograms
TH1F * m_occ_L4U
occupancy distribution for L4 U-side sensors
TH1F * m_occ_L6V
occupancy distribution for L6 V-side sensors
float m_distr_max
max of occupancy distributions
float m_distr_min
min of occupancy distributions
TH1F * m_occ_L3V
occupancy distribution for L3 V-side sensors
StoreObjPtr< SoftwareTriggerResult > m_resultStoreObjectPointer
Store Object for reading the trigger decision.
TH1F * m_occ_L3U
occupancy distribution for L3 U-side sensors
TFile * m_rootFilePtr
pointer at root file used for storing histograms
SVDHistograms< TH1F > * m_histo_zsOcc
occupancy VS ZScut histograms
StoreArray< SVDShaperDigit > m_svdShapers
SVDShaperDigit StoreArray.
TH1F * m_occ_L5V
occupancy distribution for L5 V-side sensors
std::array< APVFloatSampleType, c_nAPVSamples > APVFloatSamples
array of APVFloatSampleType objects
class to summarize SVD quantities per sensor and side
float getValue(const VxdID &vxdID, int view)
get the value contained in the corresponding bin, given VxdID and view
void fill(int layer, int ladder, int sensor, int view, float value)
fill the histogram for
TH2F * getHistogram(int view)
get a reference to the histogram for
static bool getFinalTriggerDecision(const SoftwareTriggerResult &result, bool forgetTotalResult=false)
Calculate the final cut decision using all "total_results" of all sub triggers in the software trigge...
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:39
const std::vector< VxdID > getListOfSensors() const
Get list of all sensors.
Definition: GeoCache.cc:59
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:67
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
Base class to provide Sensor Information for PXD and SVD.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
unsigned short baseType
The base integer type for VxdID.
Definition: VxdID.h:36
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
Abstract base class for different kinds of events.
STL namespace.