Belle II Software development
SVDOccupancyCalibrationsCollectorModule.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#include <svd/modules/svdOccupancyCalibrationsCollector/SVDOccupancyCalibrationsCollectorModule.h>
9
10#include <TH2F.h>
11
12using namespace std;
13using namespace Belle2;
14
15//-----------------------------------------------------------------
16// Register the Module
17//-----------------------------------------------------------------
18REG_MODULE(SVDOccupancyCalibrationsCollector);
19
20//-----------------------------------------------------------------
21// Implementation
22//-----------------------------------------------------------------
23
25{
26 //Set module properties
27
28 setDescription("This module collects hits from shaper digits to compute per sensor SVD occupancy ");
30
31 addParam("SVDShaperDigitsName", m_svdShaperDigitName, "Name of the SVDClusters list", std::string("SVDShaperDigits"));
32
33 addParam("HistogramTree", m_tree, "Name of the tree in which the histograms are saved", std::string("tree"));
34}
35
37{
38
39 m_eventMetaData.isRequired();
41
42 //create histograms
43
44 TH1F hOccupancy768("Occupancy768_L@layerL@ladderS@sensor@view", "Strip Occupancy of @layer.@ladder.@sensor @view/@side side", 768,
45 0,
46 768);
47 hOccupancy768.GetXaxis()->SetTitle("cellID");
48 TH1F hOccupancy512("Occupancy512_L@layerL@ladderS@sensor@view", "Strip Occupancy of @layer.@ladder.@sensor @view/@side side", 512,
49 0,
50 512);
51 hOccupancy512.GetXaxis()->SetTitle("cellID");
52 hm_occupancy = new SVDHistograms<TH1F>(hOccupancy768, hOccupancy768, hOccupancy768, hOccupancy512);
53
54 m_histogramTree = new TTree("tree", "tree");
55 m_histogramTree->Branch("hist", "TH1F", &m_hist, 32000, 0);
56 m_histogramTree->Branch("layer", &m_layer, "layer/I");
57 m_histogramTree->Branch("ladder", &m_ladder, "ladder/I");
58 m_histogramTree->Branch("sensor", &m_sensor, "sensor/I");
59 m_histogramTree->Branch("view", &m_side, "view/I");
60
61 m_hnevents = new TH1F("hnevents", "Number of events", 1, 0, 1);
62
63 //register objects needed to collect input to fill payloads
64 registerObject<TTree>("HTreeOccupancyCalib", m_histogramTree);
65
66}
67
69{
70
72
73 for (auto layer : geoCache.getLayers(VXD::SensorInfoBase::SVD)) {
74 for (auto ladder : geoCache.getLadders(layer)) {
75 for (Belle2::VxdID sensor : geoCache.getSensors(ladder)) {
76 for (int view = SVDHistograms<TH2F>::VIndex ; view < SVDHistograms<TH2F>::UIndex + 1; view++) {
77 // std::string s = std::string(sensor);
78 // std::string v = std::to_string(view);
79 // std::string name = string("eventT0vsCog_")+s+string("_")+v;
80 // registerObject<TH2F>(name.c_str(),hm_occupancy->getHistogram(sensor, view));
81 (hm_occupancy->getHistogram(sensor, view))->Reset();
82 }
83 }
84 }
85 }
86 m_hnevents->Reset();
87
88}
89
90
92{
93
94 int nDigits = m_storeDigits.getEntries();
95 m_hnevents->Fill(0.0); // check if HLT did not filter out the event (no rawSVD)
96
97 if (nDigits == 0)
98 return;
99
100 //loop over the SVDShaperDigits
101 int i = 0;
102 while (i < nDigits) {
103 VxdID theVxdID = m_storeDigits[i]->getSensorID();
104 int side = m_storeDigits[i]->isUStrip();
105 int CellID = m_storeDigits[i]->getCellID();
106
107 hm_occupancy->fill(theVxdID, side, CellID);
108
109 i++;
110 }
111
112
113}
114
116{
117}
118
120{
121
122 int nevents = m_hnevents->GetEntries(); //number of events processed in events
123 //getObjectPtr<TH1F>("HNevents")->GetEntries(); //number of events processed in events
124
125 B2RESULT("number of events " << nevents);
126
128 std::set<Belle2::VxdID> svdLayers = aGeometry.getLayers(VXD::SensorInfoBase::SVD);
129 std::set<Belle2::VxdID>::iterator itSvdLayers = svdLayers.begin();
130 //int itsensor = 0; //sensor numbering
131 while ((itSvdLayers != svdLayers.end())
132 && (itSvdLayers->getLayerNumber() != 7)) { //loop on Layers
133
134 std::set<Belle2::VxdID> svdLadders = aGeometry.getLadders(*itSvdLayers);
135 std::set<Belle2::VxdID>::iterator itSvdLadders = svdLadders.begin();
136
137 while (itSvdLadders != svdLadders.end()) { //loop on Ladders
138
139 std::set<Belle2::VxdID> svdSensors = aGeometry.getSensors(*itSvdLadders);
140 std::set<Belle2::VxdID>::iterator itSvdSensors = svdSensors.begin();
141
142 while (itSvdSensors != svdSensors.end()) { //loop on sensors
143
144 for (int k = 0; k < m_nSides; k ++) { //loop on Sides , k = isU(), k=0 is v-side, k=1 is u-side
145
146 (hm_occupancy->getHistogram(*itSvdSensors, k))->Scale(1. / nevents);
147 B2INFO("occupancy histo scaled by the number of events");
148 m_hist = hm_occupancy->getHistogram(*itSvdSensors, k);
149 m_layer = itSvdSensors->getLayerNumber();
150 m_ladder = itSvdSensors->getLadderNumber();
151 m_sensor = itSvdSensors->getSensorNumber();
152 m_side = k;
153
154 getObjectPtr<TTree>("HTreeOccupancyCalib")->Fill();
155 B2INFO("Filled sensors:" << m_layer << "." << m_ladder << "." << m_sensor << "." << m_side);
156 }
157 ++itSvdSensors;
158 }
159 ++itSvdLadders;
160 }
161 ++itSvdLayers;
162 }
163
164
165}
Calibration collector module base class.
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
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
StoreArray< SVDShaperDigit > m_storeDigits
shaper digits store array
TH1F * m_hnevents
first bin of the histogram is counting the processed events
TTree * m_histogramTree
tree containing as events the histograms per layer, ladder, sensor, side
void startRun() override final
Called when entering a new run.
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 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
Abstract base class for different kinds of events.
STL namespace.