Belle II Software  release-05-02-19
SVDClusterCalibrationsMonitorModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2013 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Laura Zani, Giulia Casarosa *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <svd/modules/svdCalibration/SVDClusterCalibrationsMonitorModule.h>
12 #include <vxd/geometry/GeoCache.h>
13 #include <framework/datastore/StoreObjPtr.h>
14 #include <framework/dataobjects/EventMetaData.h>
15 
16 using namespace Belle2;
17 
18 //-----------------------------------------------------------------
19 // Register the Module
20 //-----------------------------------------------------------------
21 REG_MODULE(SVDClusterCalibrationsMonitor)
22 
23 //-----------------------------------------------------------------
24 // Implementation
25 //-----------------------------------------------------------------
26 
28 {
29  // Set module properties
30  setDescription("Module to produce a list of histograms showing the uploaded calibration constants");
31 
32  // Parameter definitions
33  addParam("outputFileName", m_rootFileName, "Name of output root file.", std::string("SVDClusterCalibrationMonitor_output.root"));
34 }
35 
37 {
38 
39  m_rootFilePtr = new TFile(m_rootFileName.c_str(), "RECREATE");
40  //tree initialization
41  m_tree = new TTree("calibCls", "RECREATE");
42  b_run = m_tree->Branch("run", &m_run, "run/i");
43  b_layer = m_tree->Branch("layer", &m_layer, "layer/i");
44  b_ladder = m_tree->Branch("ladder", &m_ladder, "ladder/i");
45  b_sensor = m_tree->Branch("sensor", &m_sensor, "sensor/i");
46  b_side = m_tree->Branch("side", &m_side, "side/i");
47  b_clsSNR = m_tree->Branch("clsSNR", &m_clsSNR, "clsSNR/F");
48  b_clsSeedSNR = m_tree->Branch("clsSeedSNR", &m_clsSeedSNR, "clsSeedSNR/F");
49  b_clsAdjSNR = m_tree->Branch("clsAdjSNR", &m_clsAdjSNR, "clsAdjSNR/F");
50  b_clsScaleErr1 = m_tree->Branch("clsScaleErr1", &m_clsScaleErr1, "clsScaleErr1/F");
51  b_clsScaleErr2 = m_tree->Branch("clsScaleErr2", &m_clsScaleErr2, "clsScaleErr2/F");
52  b_clsScaleErr3 = m_tree->Branch("clsScaleErr3", &m_clsScaleErr3, "clsScaleErr3/F");
53  b_clsTimeFunc = m_tree->Branch("clsTimeFunc", &m_clsTimeFunc, "clsTimeFunc/i");
54  b_clsTimeMin = m_tree->Branch("clsTimeMin", &m_clsTimeMin, "clsTimeMin/F");
55 
56 
57  if (! m_ClusterCal.isValid())
58  B2WARNING("No valid SVDClusterCalibrations for the requested IoV");
59 
60 
62  TH1F hClsSNR("clusterMinSNR__L@layerL@ladderS@sensor@view",
63  "Cluster minimum SNR in @layer.@ladder.@sensor @view/@side",
64  100, -0.5, 99.5);
65  hClsSNR.GetXaxis()->SetTitle("cls min SNR");
66  m_hClsSNR = new SVDHistograms<TH1F>(hClsSNR);
67 
68  TH1F hClsSeedSNR("clusterSeedSNR__L@layerL@ladderS@sensor@view",
69  "Cluster Seed minimum SNR in @layer.@ladder.@sensor @view/@side",
70  100, -0.5, 99.5);
71  hClsSeedSNR.GetXaxis()->SetTitle("cls seed SNR");
72  m_hClsSeedSNR = new SVDHistograms<TH1F>(hClsSeedSNR);
73 
74 
75 
76  TH1F hClsAdjSNR("clusterAdjSNR__L@layerL@ladderS@sensor@view",
77  "Cluster Adj minimum SNR in @layer.@ladder.@sensor @view/@side",
78  100, -0.5, 99.5);
79  hClsAdjSNR.GetXaxis()->SetTitle("cls adj SNR");
80  m_hClsAdjSNR = new SVDHistograms<TH1F>(hClsAdjSNR);
81 
82 
83  //CLUSTER POSITION ERROR SCALE FACTORS
84  TH1F hClsScaleErr1("clusterScaleErr1__L@layerL@ladderS@sensor@view",
85  "Cluster Position Error Scale Factor for Size 1 in @layer.@ladder.@sensor @view/@side",
86  100, -0.5, 9.5);
87  hClsScaleErr1.GetXaxis()->SetTitle("scale factor");
88  m_hClsScaleErr1 = new SVDHistograms<TH1F>(hClsScaleErr1);
89 
90  TH1F hClsScaleErr2("clusterScaleErr2__L@layerL@ladderS@sensor@view",
91  "Cluster Position Error Scale Factor for Size 2 in @layer.@ladder.@sensor @view/@side",
92  100, -0.5, 9.5);
93  hClsScaleErr2.GetXaxis()->SetTitle("scale factor");
94  m_hClsScaleErr2 = new SVDHistograms<TH1F>(hClsScaleErr2);
95 
96  TH1F hClsScaleErr3("clusterScaleErr3__L@layerL@ladderS@sensor@view",
97  "Cluster Position Error Scale Factor for Size 3 in @layer.@ladder.@sensor @view/@side",
98  100, -0.5, 9.5);
99  hClsScaleErr3.GetXaxis()->SetTitle("scale factor");
100  m_hClsScaleErr3 = new SVDHistograms<TH1F>(hClsScaleErr3);
101 
102  //CLUSTER TIME CUTS
103  TH1F hClsTimeFuncVersion("clusterTimeSelFunction__L@layerL@ladderS@sensor@view",
104  "Cluster Time Selection Function Version in @layer.@ladder.@sensor @view/@side",
105  5, -0.5, 4.5);
106  hClsTimeFuncVersion.GetXaxis()->SetTitle("cls time selection function ID");
107  m_hClsTimeFuncVersion = new SVDHistograms<TH1F>(hClsTimeFuncVersion);
108 
109  //CLUSTER TIME CUTS
110  TH1F hClsTimeMin("clusterMinTimeS__L@layerL@ladderS@sensor@view",
111  "Cluster Minimum Time in @layer.@ladder.@sensor @view/@side",
112  200, -100.5, 99.5);
113  hClsTimeMin.GetXaxis()->SetTitle("cls min time (ns)");
114  m_hClsTimeMin = new SVDHistograms<TH1F>(hClsTimeMin);
115 
116 
117 
118 }
119 
121 {
122 
124  m_run = meta->getRun();
125 
126  if (!m_ClusterCal.isValid())
127  return;
128 
129  //call for a geometry instance
131  std::set<Belle2::VxdID> svdLayers = aGeometry.getLayers(VXD::SensorInfoBase::SVD);
132  std::set<Belle2::VxdID>::iterator itSvdLayers = svdLayers.begin();
133 
134  while ((itSvdLayers != svdLayers.end()) && (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  B2DEBUG(1, " svd sensor info " << * (svdSensors.begin()));
144 
145  while (itSvdSensors != svdSensors.end()) { //loop on sensors
146  B2DEBUG(1, " svd sensor info " << *itSvdSensors);
147 
148  m_layer = itSvdSensors->getLayerNumber();
149  m_ladder = itSvdSensors->getLadderNumber();
150  m_sensor = itSvdSensors->getSensorNumber();
152 
153  for (m_side = 0; m_side < 2; m_side++) {
154 
156  m_hClsSNR->fill(theVxdID, m_side, m_clsSNR);
157 
159  m_hClsSeedSNR->fill(theVxdID, m_side, m_clsSeedSNR);
160 
162  m_hClsAdjSNR->fill(theVxdID, m_side, m_clsAdjSNR);
163 
166 
171 
173  m_hClsTimeMin->fill(theVxdID, m_side, m_clsTimeMin);
174 
176  m_hClsTimeMin->fill(theVxdID, m_side, m_clsTimeFunc);
177 
178  m_tree->Fill();
179 
180  }
181  ++itSvdSensors;
182  }
183  ++itSvdLadders;
184  }
185  ++itSvdLayers;
186  }
187 
188 }
189 
191 {
192  B2RESULT("******************************************");
193  B2RESULT("** UNIQUE IDs of calibration DB objects **");
194  B2RESULT("");
195  if (m_ClusterCal.isValid())
196  B2RESULT(" - SVDClusterCalibrations:" << m_ClusterCal.getUniqueID());
197  else
198  B2WARNING("No valid SVDClusterCalibrations for the requested IoV");
199 
200  if (m_rootFilePtr != NULL) {
201 
202  m_rootFilePtr->cd();
203 
204  //write the tree
205  m_tree->Write();
206 
207 
208  m_rootFilePtr->mkdir("snr_cuts");
209  m_rootFilePtr->mkdir("time_cuts");
210  m_rootFilePtr->mkdir("scale_factor");
211 
213 
214  for (auto layer : geoCache.getLayers(VXD::SensorInfoBase::SVD))
215  for (auto ladder : geoCache.getLadders(layer))
216  for (Belle2::VxdID sensor : geoCache.getSensors(ladder))
217  for (int view = SVDHistograms<TH1F>::VIndex ; view < SVDHistograms<TH1F>::UIndex + 1; view++) {
218 
219  //writing the histograms to root:
220 
221  m_rootFilePtr->cd("snr_cuts");
222  (m_hClsSNR->getHistogram(sensor, view))->Write();
223  (m_hClsSeedSNR->getHistogram(sensor, view))->Write();
224  (m_hClsAdjSNR->getHistogram(sensor, view))->Write();
225 
226 
227  m_rootFilePtr->cd("time_cuts");
228  (m_hClsTimeFuncVersion->getHistogram(sensor, view))->Write();
229  (m_hClsTimeMin->getHistogram(sensor, view))->Write();
230 
231  m_rootFilePtr->cd("scale_factor");
232  (m_hClsScaleErr1->getHistogram(sensor, view))->Write();
233  (m_hClsScaleErr2->getHistogram(sensor, view))->Write();
234  (m_hClsScaleErr3->getHistogram(sensor, view))->Write();
235 
236 
237  }
238 
239  m_rootFilePtr->Close();
240  B2RESULT("The rootfile containing the list of histograms has been filled and closed [Cluster].");
241 
242 
243  }
244 }
245 
Belle2::SVDClusterCalibrationsMonitorModule::m_ladder
int m_ladder
ladder number
Definition: SVDClusterCalibrationsMonitorModule.h:89
Belle2::VxdID
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:43
Belle2::SVDClusterCalibrationsMonitorModule::m_hClsScaleErr2
SVDHistograms< TH1F > * m_hClsScaleErr2
cluster size 2, position error scale factor histo
Definition: SVDClusterCalibrationsMonitorModule.h:112
Belle2::SVDClusterCalibrationsMonitorModule::endRun
virtual void endRun() override
print the payloads uniqueID and write tree and histograms to the rootfile
Definition: SVDClusterCalibrationsMonitorModule.cc:190
Belle2::SVDClusterCalibrationsMonitorModule::m_hClsSeedSNR
SVDHistograms< TH1F > * m_hClsSeedSNR
cluster seed SNR histo
Definition: SVDClusterCalibrationsMonitorModule.h:107
Belle2::SVDClusterCalibrationsMonitorModule::m_clsScaleErr1
float m_clsScaleErr1
cluster size 1 pos err scale factor SNR
Definition: SVDClusterCalibrationsMonitorModule.h:95
Belle2::SVDHistograms< TH1F >
Belle2::SVDClusterCalibrationsMonitorModule::b_clsTimeMin
TBranch * b_clsTimeMin
cluster cut min time
Definition: SVDClusterCalibrationsMonitorModule.h:83
Belle2::SVDClusterCalibrations::getTimeSelectionFunction
int getTimeSelectionFunction(const Belle2::VxdID &sensorID, const bool &isU) const
Return the version of the function used to determine whether the cluster time is acceptable at the SP...
Definition: SVDClusterCalibrations.h:226
Belle2::SVDClusterCalibrations::getMinSeedSNR
double getMinSeedSNR(const Belle2::VxdID &sensorID, const bool &isU) const
Return the minimum SNR for the seed.
Definition: SVDClusterCalibrations.h:99
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::SVDClusterCalibrations::getMinClusterSNR
double getMinClusterSNR(const Belle2::VxdID &sensorID, const bool &isU) const
Return the minimum SNR for the cluster.
Definition: SVDClusterCalibrations.h:146
Belle2::SVDClusterCalibrationsMonitorModule::m_rootFileName
std::string m_rootFileName
root file name
Definition: SVDClusterCalibrationsMonitorModule.h:64
Belle2::SVDClusterCalibrationsMonitorModule::b_clsScaleErr1
TBranch * b_clsScaleErr1
cluster size 1 pos err scale factor SNR
Definition: SVDClusterCalibrationsMonitorModule.h:79
Belle2::SVDClusterCalibrationsMonitorModule::b_clsAdjSNR
TBranch * b_clsAdjSNR
cluster adj SNR
Definition: SVDClusterCalibrationsMonitorModule.h:78
Belle2::SVDClusterCalibrationsMonitorModule::b_clsScaleErr3
TBranch * b_clsScaleErr3
cluster size 3 pos err scale factor SNR
Definition: SVDClusterCalibrationsMonitorModule.h:81
Belle2::SVDClusterCalibrationsMonitorModule::b_layer
TBranch * b_layer
layer number
Definition: SVDClusterCalibrationsMonitorModule.h:73
Belle2::SVDClusterCalibrationsMonitorModule::event
virtual void event() override
fill trees and histograms
Definition: SVDClusterCalibrationsMonitorModule.cc:120
Belle2::SVDClusterCalibrations::getUniqueID
TString getUniqueID()
returns the unique ID of the payload
Definition: SVDClusterCalibrations.h:267
Belle2::SVDClusterCalibrationsMonitorModule::m_sensor
int m_sensor
sensor number
Definition: SVDClusterCalibrationsMonitorModule.h:90
Belle2::SVDClusterCalibrationsMonitorModule::m_clsSNR
float m_clsSNR
cluster SNR
Definition: SVDClusterCalibrationsMonitorModule.h:92
Belle2::SVDClusterCalibrationsMonitorModule::m_hClsScaleErr1
SVDHistograms< TH1F > * m_hClsScaleErr1
cluster size 1, position error scale factor histo
Definition: SVDClusterCalibrationsMonitorModule.h:111
Belle2::VXD::GeoCache::getLayers
const std::set< Belle2::VxdID > getLayers(SensorInfoBase::SensorType sensortype=SensorInfoBase::VXD)
Return a set of all known Layers.
Definition: GeoCache.cc:177
Belle2::SVDClusterCalibrations::getMinClusterTime
float getMinClusterTime(const Belle2::VxdID &sensorID, const bool &isU) const
Return the min value of the cluster time to use it for reconstruction.
Definition: SVDClusterCalibrations.h:251
Belle2::SVDClusterCalibrationsMonitorModule::beginRun
virtual void beginRun() override
initialize the TTrees and check validities of payloads
Definition: SVDClusterCalibrationsMonitorModule.cc:36
Belle2::SVDClusterCalibrationsMonitorModule::b_clsSeedSNR
TBranch * b_clsSeedSNR
cluster seed SNR
Definition: SVDClusterCalibrationsMonitorModule.h:77
Belle2::VXD::GeoCache::getSensors
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:205
Belle2::SVDClusterCalibrationsMonitorModule::b_clsSNR
TBranch * b_clsSNR
cluster SNR
Definition: SVDClusterCalibrationsMonitorModule.h:76
Belle2::Module
Base class for Modules.
Definition: Module.h:74
Belle2::SVDClusterCalibrationsMonitorModule
Module to produce a list of histogram showing the uploaded calibration constants.
Definition: SVDClusterCalibrationsMonitorModule.h:45
Belle2::SVDClusterCalibrationsMonitorModule::m_ClusterCal
SVDClusterCalibrations m_ClusterCal
cluster calibrations payload
Definition: SVDClusterCalibrationsMonitorModule.h:103
Belle2::SVDClusterCalibrationsMonitorModule::b_ladder
TBranch * b_ladder
ladder number
Definition: SVDClusterCalibrationsMonitorModule.h:72
Belle2::SVDClusterCalibrationsMonitorModule::m_side
int m_side
sensor side
Definition: SVDClusterCalibrationsMonitorModule.h:91
Belle2::VXD::GeoCache::getInstance
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:215
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::StoreObjPtr
Type-safe access to single objects in the data store.
Definition: ParticleList.h:33
Belle2::SVDClusterCalibrationsMonitorModule::m_hClsAdjSNR
SVDHistograms< TH1F > * m_hClsAdjSNR
cluster adj SNR histo
Definition: SVDClusterCalibrationsMonitorModule.h:108
Belle2::SVDClusterCalibrations::isValid
bool isValid()
returns true if the m_aDBObtPtr is valid in the requested IoV
Definition: SVDClusterCalibrations.h:270
Belle2::SVDClusterCalibrationsMonitorModule::m_clsTimeFunc
int m_clsTimeFunc
cluster cut time function version
Definition: SVDClusterCalibrationsMonitorModule.h:98
Belle2::SVDClusterCalibrationsMonitorModule::b_sensor
TBranch * b_sensor
sensor number
Definition: SVDClusterCalibrationsMonitorModule.h:74
Belle2::SVDClusterCalibrations::getMinAdjSNR
double getMinAdjSNR(const Belle2::VxdID &sensorID, const bool &isU) const
Return the minimum SNR for the adjacent.
Definition: SVDClusterCalibrations.h:122
Belle2::SVDClusterCalibrationsMonitorModule::b_run
TBranch * b_run
run number
Definition: SVDClusterCalibrationsMonitorModule.h:71
Belle2::VXD::SensorInfoBase::SVD
@ SVD
SVD Sensor.
Definition: SensorInfoBase.h:45
Belle2::SVDClusterCalibrationsMonitorModule::m_hClsTimeMin
SVDHistograms< TH1F > * m_hClsTimeMin
cluster cut minimum time histo
Definition: SVDClusterCalibrationsMonitorModule.h:117
Belle2::SVDClusterCalibrationsMonitorModule::m_clsAdjSNR
float m_clsAdjSNR
cluster adj SNR
Definition: SVDClusterCalibrationsMonitorModule.h:94
Belle2::SVDHistograms::fill
void fill(const VxdID &vxdID, int view, Types ... args)
fill the histogram for
Definition: SVDHistograms.h:89
Belle2::SVDClusterCalibrationsMonitorModule::m_run
int m_run
run number
Definition: SVDClusterCalibrationsMonitorModule.h:87
Belle2::VXD::GeoCache
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:41
Belle2::SVDClusterCalibrationsMonitorModule::m_clsScaleErr3
float m_clsScaleErr3
cluster size 3 pos err scale factor SNR
Definition: SVDClusterCalibrationsMonitorModule.h:97
Belle2::SVDClusterCalibrationsMonitorModule::m_rootFilePtr
TFile * m_rootFilePtr
pointer at root file used for storing histograms
Definition: SVDClusterCalibrationsMonitorModule.h:67
Belle2::SVDHistograms::getHistogram
H * getHistogram(const VxdID &vxdID, int view)
get a reference to the histogram for
Definition: SVDHistograms.h:68
Belle2::SVDClusterCalibrationsMonitorModule::b_clsScaleErr2
TBranch * b_clsScaleErr2
cluster size 2 pos err scale factor SNR
Definition: SVDClusterCalibrationsMonitorModule.h:80
Belle2::SVDClusterCalibrationsMonitorModule::b_clsTimeFunc
TBranch * b_clsTimeFunc
cluster cut time function ID
Definition: SVDClusterCalibrationsMonitorModule.h:82
Belle2::SVDClusterCalibrationsMonitorModule::m_hClsScaleErr3
SVDHistograms< TH1F > * m_hClsScaleErr3
cluster size 3, position error scale factor histo
Definition: SVDClusterCalibrationsMonitorModule.h:113
Belle2::VXD::GeoCache::getLadders
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:194
Belle2::SVDClusterCalibrationsMonitorModule::m_hClsTimeFuncVersion
SVDHistograms< TH1F > * m_hClsTimeFuncVersion
cluster cut time function version histo
Definition: SVDClusterCalibrationsMonitorModule.h:116
Belle2::SVDClusterCalibrationsMonitorModule::m_layer
int m_layer
layer number
Definition: SVDClusterCalibrationsMonitorModule.h:88
Belle2::SVDClusterCalibrationsMonitorModule::b_side
TBranch * b_side
sensor side
Definition: SVDClusterCalibrationsMonitorModule.h:75
Belle2::SVDClusterCalibrationsMonitorModule::m_hClsSNR
SVDHistograms< TH1F > * m_hClsSNR
cluster SNR histo
Definition: SVDClusterCalibrationsMonitorModule.h:106
Belle2::SVDClusterCalibrations::getCorrectedClusterPositionError
double getCorrectedClusterPositionError(const Belle2::VxdID &sensorID, const bool &isU, const int &size, const double &raw_error) const
Return the corrected cluster position error.
Definition: SVDClusterCalibrations.h:75
Belle2::SVDClusterCalibrationsMonitorModule::m_tree
TTree * m_tree
pointer at tree containing the mean and RMS of calibration constants
Definition: SVDClusterCalibrationsMonitorModule.h:68
Belle2::SVDClusterCalibrationsMonitorModule::m_clsSeedSNR
float m_clsSeedSNR
cluster seed SNR
Definition: SVDClusterCalibrationsMonitorModule.h:93
Belle2::SVDClusterCalibrationsMonitorModule::m_clsScaleErr2
float m_clsScaleErr2
cluster size 2 pos err scale factor SNR
Definition: SVDClusterCalibrationsMonitorModule.h:96
Belle2::SVDClusterCalibrationsMonitorModule::m_clsTimeMin
float m_clsTimeMin
cluster cut min time
Definition: SVDClusterCalibrationsMonitorModule.h:99