8 #include <svd/modules/svdOccupancyCalibrationsCollector/SVDOccupancyCalibrationsCollectorModule.h>
28 setDescription(
"This module collects hits from shaper digits to compute per sensor SVD occupancy ");
29 setPropertyFlags(c_ParallelProcessingCertified);
31 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"));
36 void SVDOccupancyCalibrationsCollectorModule::prepare()
39 m_eventMetaData.isRequired();
40 m_storeDigits.isRequired(m_svdShaperDigitName);
44 TH1F hOccupancy768(
"Occupancy768_L@layerL@ladderS@sensor@view",
"Strip Occupancy of @layer.@ladder.@sensor @view/@side side", 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,
51 hOccupancy512.GetXaxis()->SetTitle(
"cellID");
52 hm_occupancy =
new SVDHistograms<TH1F>(hOccupancy768, hOccupancy768, hOccupancy768, hOccupancy512);
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");
61 m_hnevents =
new TH1F(
"hnevents",
"Number of events", 1, 0, 1);
64 registerObject<TTree>(
"HTreeOccupancyCalib", m_histogramTree);
68 void SVDOccupancyCalibrationsCollectorModule::startRun()
73 for (
auto layer : geoCache.getLayers(VXD::SensorInfoBase::SVD)) {
74 for (
auto ladder : geoCache.getLadders(layer)) {
81 (hm_occupancy->getHistogram(sensor, view))->Reset();
91 void SVDOccupancyCalibrationsCollectorModule::collect()
94 int nDigits = m_storeDigits.getEntries();
95 m_hnevents->Fill(0.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();
107 hm_occupancy->fill(theVxdID, side, CellID);
115 void SVDOccupancyCalibrationsCollectorModule::finish()
119 void SVDOccupancyCalibrationsCollectorModule::closeRun()
122 int nevents = m_hnevents->GetEntries();
125 B2RESULT(
"number of events " << nevents);
128 std::set<Belle2::VxdID> svdLayers = aGeometry.
getLayers(VXD::SensorInfoBase::SVD);
129 std::set<Belle2::VxdID>::iterator itSvdLayers = svdLayers.begin();
131 while ((itSvdLayers != svdLayers.end())
132 && (itSvdLayers->getLayerNumber() != 7)) {
134 std::set<Belle2::VxdID> svdLadders = aGeometry.
getLadders(*itSvdLayers);
135 std::set<Belle2::VxdID>::iterator itSvdLadders = svdLadders.begin();
137 while (itSvdLadders != svdLadders.end()) {
139 std::set<Belle2::VxdID> svdSensors = aGeometry.
getSensors(*itSvdLadders);
140 std::set<Belle2::VxdID>::iterator itSvdSensors = svdSensors.begin();
142 while (itSvdSensors != svdSensors.end()) {
144 for (
int k = 0; k < m_nSides; k ++) {
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();
154 getObjectPtr<TTree>(
"HTreeOccupancyCalib")->Fill();
155 B2INFO(
"Filled sensors:" << m_layer <<
"." << m_ladder <<
"." << m_sensor <<
"." << m_side);
Calibration collector module base class.
This This module collects hits from shaper digits to compute per sensor SVD occupancy using mu+mu- ev...
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
const std::set< Belle2::VxdID > getLayers(SensorInfoBase::SensorType sensortype=SensorInfoBase::VXD)
Return a set of all known Layers.
const std::set< Belle2::VxdID > & getSensors(Belle2::VxdID ladder) const
Return a set of all sensor IDs belonging to a given ladder.
const std::set< Belle2::VxdID > & getLadders(Belle2::VxdID layer) const
Return a set of all ladder IDs belonging to a given layer.
Class to uniquely identify a any structure of the PXD and SVD.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Abstract base class for different kinds of events.