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#include <hlt/softwaretrigger/core/FinalTriggerDecisionCalculator.h>
10
11#include <TH2F.h>
12
13using namespace std;
14using namespace Belle2;
15
16//-----------------------------------------------------------------
17// Register the Module
18//-----------------------------------------------------------------
19REG_MODULE(SVDOccupancyCalibrationsCollector);
20
21//-----------------------------------------------------------------
22// Implementation
23//-----------------------------------------------------------------
24
26{
27 //Set module properties
28
29 setDescription("This module collects hits from shaper digits to compute per sensor SVD occupancy ");
31
32 addParam("SVDShaperDigitsName", m_svdShaperDigitName, "Name of the SVDClusters list", std::string("SVDShaperDigits"));
33 addParam("HistogramTree", m_tree, "Name of the tree in which the histograms are saved", std::string("tree"));
34 addParam("skipHLTRejectedEvents", m_skipRejectedEvents, "If True, skip events rejected by HLT.", bool(false));
35}
36
38{
39
40 m_eventMetaData.isRequired();
42
43 //create histograms
44
45 TH1F hOccupancy768("Occupancy768_L@layerL@ladderS@sensor@view", "Strip Occupancy of @layer.@ladder.@sensor @view/@side side", 768,
46 0,
47 768);
48 hOccupancy768.GetXaxis()->SetTitle("cellID");
49 TH1F hOccupancy512("Occupancy512_L@layerL@ladderS@sensor@view", "Strip Occupancy of @layer.@ladder.@sensor @view/@side side", 512,
50 0,
51 512);
52 hOccupancy512.GetXaxis()->SetTitle("cellID");
53 hm_occupancy = new SVDHistograms<TH1F>(hOccupancy768, hOccupancy768, hOccupancy768, hOccupancy512);
54
55 m_histogramTree = new TTree("tree", "tree");
56 m_histogramTree->Branch("hist", "TH1F", &m_hist, 32000, 0);
57 m_histogramTree->Branch("layer", &m_layer, "layer/I");
58 m_histogramTree->Branch("ladder", &m_ladder, "ladder/I");
59 m_histogramTree->Branch("sensor", &m_sensor, "sensor/I");
60 m_histogramTree->Branch("view", &m_side, "view/I");
61
62 m_hnevents = new TH1F("hnevents", "Number of events", 3, 0, 2);
63
64 //register objects needed to collect input to fill payloads
65 registerObject<TTree>("HTreeOccupancyCalib", m_histogramTree);
66 registerObject<TH1F>("HNEvents", m_hnevents);
67
68}
69
71{
72
74
75 for (auto layer : geoCache.getLayers(VXD::SensorInfoBase::SVD)) {
76 for (auto ladder : geoCache.getLadders(layer)) {
77 for (Belle2::VxdID sensor : geoCache.getSensors(ladder)) {
78 for (int view = SVDHistograms<TH2F>::VIndex ; view < SVDHistograms<TH2F>::UIndex + 1; view++) {
79 (hm_occupancy->getHistogram(sensor, view))->Reset();
80 }
81 }
82 }
83 }
84 // m_hnevents->Reset();
85 getObjectPtr<TH1F>("HNEvents")->Reset();
86
87}
88
89
91{
92
95 if (!eventAccepted) return;
96 }
97
98 int nDigits = m_storeDigits.getEntries();
99 getObjectPtr<TH1F>("HNEvents")->Fill(1);
100
101 if (nDigits == 0)
102 return;
103
104 //loop over the SVDShaperDigits
105 int i = 0;
106 while (i < nDigits) {
107 VxdID theVxdID = m_storeDigits[i]->getSensorID();
108 int side = m_storeDigits[i]->isUStrip();
109 int CellID = m_storeDigits[i]->getCellID();
110
111 hm_occupancy->fill(theVxdID, side, CellID);
112
113 i++;
114 }
115
116
117}
118
120{
121}
122
124{
125
126 int nevents = getObjectPtr<TH1F>("HNEvents")->GetEntries(); //number of events processed in events
127
128 B2RESULT("number of events " << nevents);
129
131 std::set<Belle2::VxdID> svdLayers = aGeometry.getLayers(VXD::SensorInfoBase::SVD);
132 std::set<Belle2::VxdID>::iterator itSvdLayers = svdLayers.begin();
133 while ((itSvdLayers != svdLayers.end())
134 && (itSvdLayers->getLayerNumber() != 7)) { //loop on Layers
135
136 std::set<Belle2::VxdID> svdLadders = aGeometry.getLadders(*itSvdLayers);
137 std::set<Belle2::VxdID>::iterator itSvdLadders = svdLadders.begin();
138
139 while (itSvdLadders != svdLadders.end()) { //loop on Ladders
140
141 std::set<Belle2::VxdID> svdSensors = aGeometry.getSensors(*itSvdLadders);
142 std::set<Belle2::VxdID>::iterator itSvdSensors = svdSensors.begin();
143
144 while (itSvdSensors != svdSensors.end()) { //loop on sensors
145
146 for (int k = 0; k < m_nSides; k ++) { //loop on Sides , k = isU(), k=0 is v-side, k=1 is u-side
147
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
bool m_skipRejectedEvents
if true skip events rejected by HLT (default)
TH1F * m_hnevents
first bin of the histogram is counting the processed events
StoreObjPtr< SoftwareTriggerResult > m_resultStoreObjectPointer
Store Object for reading the trigger decision.
TTree * m_histogramTree
tree containing as events the histograms per layer, ladder, sensor, side
void startRun() override final
Called when entering a new run.
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 facilitate easy access to sensor information of the VXD like coordinate transformations or p...
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:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:649
Abstract base class for different kinds of events.
STL namespace.