Belle II Software development
SVDTimeCalibrationsMonitorModule.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/svdCalibration/SVDTimeCalibrationsMonitorModule.h>
10#include <vxd/geometry/GeoCache.h>
11#include <framework/datastore/StoreObjPtr.h>
12#include <framework/dataobjects/EventMetaData.h>
13
14using namespace Belle2;
15
16//-----------------------------------------------------------------
17// Register the Module
18//-----------------------------------------------------------------
19REG_MODULE(SVDTimeCalibrationsMonitor);
20
21//-----------------------------------------------------------------
22// Implementation
23//-----------------------------------------------------------------
24
26{
27 // Set module properties
28 setDescription("Module to produce a list of histograms showing the uploaded calibration constants");
29
30 // Parameter definitions
31 addParam("outputFileName", m_rootFileName, "Name of output root file.", std::string("SVDTimeCalibrationMonitor_output.root"));
32
33 addParam("timeAlgo", m_timeAlgo, "Name of the time algorithm: CoG6, CoG3, ELS3.", std::string("CoG3"));
34}
35
37{
38
39 m_rootFilePtr = new TFile(m_rootFileName.c_str(), "RECREATE");
40 //tree initialization
41 m_tree = new TTree("calibTime", "RECREATE");
42 b_exp = m_tree->Branch("exp", &m_exp, "exp/i");
43 b_run = m_tree->Branch("run", &m_run, "run/i");
44 b_layer = m_tree->Branch("layer", &m_layer, "layer/i");
45 b_ladder = m_tree->Branch("ladder", &m_ladder, "ladder/i");
46 b_sensor = m_tree->Branch("sensor", &m_sensor, "sensor/i");
47 b_side = m_tree->Branch("side", &m_side, "side/i");
48 b_triggerBin = m_tree->Branch("triggerBin", &m_triggerBin, "triggerBin/i");
49 b_c0 = m_tree->Branch("c0", &m_c0, "c0/F");
50 b_c1 = m_tree->Branch("c1", &m_c1, "c1/F");
51 b_c2 = m_tree->Branch("c2", &m_c2, "c2/F");
52 b_c3 = m_tree->Branch("c3", &m_c3, "c3/F");
53
54
55 if (TString(m_timeAlgo).Contains("CoG6") && ! m_CoG6TimeCal.isValid())
56 B2WARNING("No valid SVDCoGTimeCalibrations for the requested IoV");
57
58 if (TString(m_timeAlgo).Contains("CoG3") && ! m_CoG3TimeCal.isValid())
59 B2WARNING("No valid SVD3SampleCoGTimeCalibrations for the requested IoV");
60
61
62 if (TString(m_timeAlgo).Contains("ELS3") && ! m_ELS3TimeCal.isValid())
63 B2WARNING("No valid SVD3SampleELSTimeCalibrations for the requested IoV");
64
65
66}
67
69{
70
72 m_run = meta->getRun();
73 m_exp = meta->getExperiment();
74
75 if (TString(m_timeAlgo).Contains("CoG6"))
77 return;
78
79 if (TString(m_timeAlgo).Contains("CoG3"))
81 return;
82
83 if (TString(m_timeAlgo).Contains("ELS3"))
85 return;
86
87
88 //call for a geometry instance
90 std::set<Belle2::VxdID> svdLayers = aGeometry.getLayers(VXD::SensorInfoBase::SVD);
91 std::set<Belle2::VxdID>::iterator itSvdLayers = svdLayers.begin();
92
93 while ((itSvdLayers != svdLayers.end()) && (itSvdLayers->getLayerNumber() != 7)) { //loop on Layers
94
95 std::set<Belle2::VxdID> svdLadders = aGeometry.getLadders(*itSvdLayers);
96 std::set<Belle2::VxdID>::iterator itSvdLadders = svdLadders.begin();
97
98 while (itSvdLadders != svdLadders.end()) { //loop on Ladders
99
100 std::set<Belle2::VxdID> svdSensors = aGeometry.getSensors(*itSvdLadders);
101 std::set<Belle2::VxdID>::iterator itSvdSensors = svdSensors.begin();
102 B2DEBUG(1, " svd sensor info " << * (svdSensors.begin()));
103
104 while (itSvdSensors != svdSensors.end()) { //loop on sensors
105 B2DEBUG(1, " svd sensor info " << *itSvdSensors);
106
107 m_layer = itSvdSensors->getLayerNumber();
108 m_ladder = itSvdSensors->getLadderNumber();
109 m_sensor = itSvdSensors->getSensorNumber();
111
112 for (m_side = 0; m_side < 2; m_side++) {
113
114 for (m_triggerBin = 0; m_triggerBin < 4; m_triggerBin++) {
115
116 if (TString(m_timeAlgo).Contains("CoG6")) {
117 m_CoG6TimeCal.getCorrectedTime(theVxdID, m_side, 0 /*strip*/, 0 /*raw time*/, m_triggerBin);
118 m_CoG6TimeCal.getCorrectedTime(theVxdID, m_side, 0 /*strip*/, 1 /*raw time*/, m_triggerBin);
119 m_CoG6TimeCal.getCorrectedTime(theVxdID, m_side, 0 /*strip*/, 2 /*raw time*/, m_triggerBin);
120 m_CoG6TimeCal.getCorrectedTime(theVxdID, m_side, 0 /*strip*/, 4 /*raw time*/, m_triggerBin);
121 }
122
123 if (TString(m_timeAlgo).Contains("CoG3")) {
124 double f0 = m_CoG3TimeCal.getCorrectedTime(theVxdID, m_side, 0 /*strip*/, 0 /*raw time*/, m_triggerBin);
125 double f1 = m_CoG3TimeCal.getCorrectedTime(theVxdID, m_side, 0 /*strip*/, 1 /*raw time*/, m_triggerBin);
126 double f2 = m_CoG3TimeCal.getCorrectedTime(theVxdID, m_side, 0 /*strip*/, 2 /*raw time*/, m_triggerBin);
127 double f4 = m_CoG3TimeCal.getCorrectedTime(theVxdID, m_side, 0 /*strip*/, 4 /*raw time*/, m_triggerBin);
128
129 m_c0 = f0;
130 m_c1 = 1. / 12 * (-21 * f0 + 32 * f1 - 12 * f2 + f4);
131 m_c2 = 1. / 8 * (7 * f0 - 16 * f1 + 10 * f2 - f4);
132 m_c3 = 1. / 24 * (-3 * f0 + 8 * f1 - 6 * f2 + f4);
133
134 }
135
136 if (TString(m_timeAlgo).Contains("ELS3")) {
137 m_ELS3TimeCal.getCorrectedTime(theVxdID, m_side, 0 /*strip*/, 0 /*raw time*/, m_triggerBin);
138 m_ELS3TimeCal.getCorrectedTime(theVxdID, m_side, 0 /*strip*/, 1 /*raw time*/, m_triggerBin);
139 m_ELS3TimeCal.getCorrectedTime(theVxdID, m_side, 0 /*strip*/, 2 /*raw time*/, m_triggerBin);
140 m_ELS3TimeCal.getCorrectedTime(theVxdID, m_side, 0 /*strip*/, 4 /*raw time*/, m_triggerBin);
141 }
142
143 m_tree->Fill();
144 }
145
146 }
147 ++itSvdSensors;
148 }
149 ++itSvdLadders;
150 }
151 ++itSvdLayers;
152 }
153
154}
155
157{
158 B2RESULT("******************************************");
159 B2RESULT("** UNIQUE IDs of calibration DB objects **");
160 B2RESULT("");
161 if (TString(m_timeAlgo).Contains("CoG6")) {
163 B2RESULT(" - SVDCoGTimeCalibrations:" << m_CoG6TimeCal.getUniqueID());
164 else
165 B2WARNING("No valid SVDCoGTimeCalibrations for the requested IoV");
166 }
167
168 if (TString(m_timeAlgo).Contains("CoG3")) {
170 B2RESULT(" - SVD3SampleCoGTimeCalibrations:" << m_CoG3TimeCal.getUniqueID());
171 else
172 B2WARNING("No valid SVD3SampleCoGTimeCalibrations for the requested IoV");
173 }
174
175 if (TString(m_timeAlgo).Contains("ELS3")) {
177 B2RESULT(" - SVD3SampleELSTimeCalibrations:" << m_ELS3TimeCal.getUniqueID());
178 else
179 B2WARNING("No valid SVD3SampleELSTimeCalibrations for the requested IoV");
180 }
181
182 if (m_rootFilePtr != nullptr) {
183
184 m_rootFilePtr->cd();
185
186 //write the tree
187 m_tree->Write();
188
189 m_rootFilePtr->Close();
190 B2RESULT("The rootfile containing the list of histograms has been filled and closed [Time].");
191
192
193 }
194}
195
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
TString getUniqueID()
returns the unique ID of the payload
double getCorrectedTime(const Belle2::VxdID &sensorID, const bool &isU, const unsigned short &strip, const double &raw_time, const int &bin) const
Return the charge (number of electrons/holes) collected on a specific strip, given the number of ADC ...
bool isValid()
returns true if the m_aDBObtPtr is valid in the requested IoV
TString getUniqueID()
returns the unique ID of the payload
double getCorrectedTime(const Belle2::VxdID &sensorID, const bool &isU, const unsigned short &strip, const double &raw_time, const int &bin) const
Return the charge (number of electrons/holes) collected on a specific strip, given the number of ADC ...
bool isValid()
returns true if the m_aDBObtPtr is valid in the requested IoV
TString getUniqueID()
returns the unique ID of the payload
double getCorrectedTime(const Belle2::VxdID &sensorID, const bool &isU, const unsigned short &strip, const double &raw_time, const int &bin) const
Return the strip time, given the raw strip time.
bool isValid()
returns true if the m_aDBObtPtr is valid in the requested IoV
TTree * m_tree
pointer at tree containing the mean and RMS of calibration constants
SVDTimeCalibrationsMonitorModule()
Constructor: Sets the description, the properties and the parameters of the module.
virtual void event() override
fill trees and histograms
virtual void endRun() override
print the payloads uniqueID and write tree and histograms to the rootfile
SVD3SampleCoGTimeCalibrations m_CoG3TimeCal
CoG3 time calibrations payload.
SVDCoGTimeCalibrations m_CoG6TimeCal
CoG6 time calibrations payload.
virtual void beginRun() override
initialize the TTrees and check validities of payloads
SVD3SampleELSTimeCalibrations m_ELS3TimeCal
ELS3 time calibrations payload.
TFile * m_rootFilePtr
pointer at root file used for storing histograms
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
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.