Belle II Software  release-05-01-25
SVDCoGTimeCalibrationsMonitorModule.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/SVDCoGTimeCalibrationsMonitorModule.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(SVDCoGTimeCalibrationsMonitor)
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("SVDCoGTimeCalibrationMonitor_output.root"));
34 }
35 
37 {
38 
39  m_rootFilePtr = new TFile(m_rootFileName.c_str(), "RECREATE");
40  //tree initialization
41  m_tree = new TTree("calibCoGTime", "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_triggerBin = m_tree->Branch("triggerBin", &m_triggerBin, "triggerBin/i");
48  b_alpha = m_tree->Branch("alpha", &m_alpha, "alpha/F");
49  b_beta = m_tree->Branch("beta", &m_beta, "beta/F");
50 
51  if (! m_CoGTimeCal.isValid())
52  B2WARNING("No valid SVDCoGTimeCalibrations for the requested IoV");
53 
54 
56 
57 
58 }
59 
61 {
62 
64  m_run = meta->getRun();
65 
66  if (!m_CoGTimeCal.isValid())
67  return;
68 
69  //call for a geometry instance
71  std::set<Belle2::VxdID> svdLayers = aGeometry.getLayers(VXD::SensorInfoBase::SVD);
72  std::set<Belle2::VxdID>::iterator itSvdLayers = svdLayers.begin();
73 
74  while ((itSvdLayers != svdLayers.end()) && (itSvdLayers->getLayerNumber() != 7)) { //loop on Layers
75 
76  std::set<Belle2::VxdID> svdLadders = aGeometry.getLadders(*itSvdLayers);
77  std::set<Belle2::VxdID>::iterator itSvdLadders = svdLadders.begin();
78 
79  while (itSvdLadders != svdLadders.end()) { //loop on Ladders
80 
81  std::set<Belle2::VxdID> svdSensors = aGeometry.getSensors(*itSvdLadders);
82  std::set<Belle2::VxdID>::iterator itSvdSensors = svdSensors.begin();
83  B2DEBUG(1, " svd sensor info " << * (svdSensors.begin()));
84 
85  while (itSvdSensors != svdSensors.end()) { //loop on sensors
86  B2DEBUG(1, " svd sensor info " << *itSvdSensors);
87 
88  m_layer = itSvdSensors->getLayerNumber();
89  m_ladder = itSvdSensors->getLadderNumber();
90  m_sensor = itSvdSensors->getSensorNumber();
92 
93  for (m_side = 0; m_side < 2; m_side++) {
94 
95  for (m_triggerBin = 0; m_triggerBin < 4; m_triggerBin++) {
96 
97  m_beta = m_CoGTimeCal.getCorrectedTime(theVxdID, m_side, 0 /*strip*/, 0 /*raw time*/, m_triggerBin);
98  m_alpha = m_CoGTimeCal.getCorrectedTime(theVxdID, m_side, 0 /*strip*/, 1 /*raw time*/, m_triggerBin) - m_beta;
99 
100  m_tree->Fill();
101  }
102 
103  }
104  ++itSvdSensors;
105  }
106  ++itSvdLadders;
107  }
108  ++itSvdLayers;
109  }
110 
111 }
112 
114 {
115  B2RESULT("******************************************");
116  B2RESULT("** UNIQUE IDs of calibration DB objects **");
117  B2RESULT("");
118  if (m_CoGTimeCal.isValid())
119  B2RESULT(" - SVDCoGTimeCalibrations:" << m_CoGTimeCal.getUniqueID());
120  else
121  B2WARNING("No valid SVDCoGTimeCalibrations for the requested IoV");
122 
123  if (m_rootFilePtr != NULL) {
124 
125  m_rootFilePtr->cd();
126 
127  //write the tree
128  m_tree->Write();
129 
130  m_rootFilePtr->Close();
131  B2RESULT("The rootfile containing the list of histograms has been filled and closed [CoGTime].");
132 
133 
134  }
135 }
136 
Belle2::SVDCoGTimeCalibrationsMonitorModule::event
virtual void event() override
fill trees and histograms
Definition: SVDCoGTimeCalibrationsMonitorModule.cc:60
Belle2::SVDCoGTimeCalibrationsMonitorModule::b_ladder
TBranch * b_ladder
ladder number
Definition: SVDCoGTimeCalibrationsMonitorModule.h:69
Belle2::VxdID
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:43
Belle2::SVDCoGTimeCalibrationsMonitorModule::m_rootFileName
std::string m_rootFileName
root file name
Definition: SVDCoGTimeCalibrationsMonitorModule.h:61
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::SVDCoGTimeCalibrations::getUniqueID
TString getUniqueID()
returns the unique ID of the payload
Definition: SVDCoGTimeCalibrations.h:112
Belle2::SVDCoGTimeCalibrationsMonitorModule
Module to produce a list of histogram showing the uploaded calibration constants.
Definition: SVDCoGTimeCalibrationsMonitorModule.h:42
Belle2::SVDCoGTimeCalibrationsMonitorModule::b_layer
TBranch * b_layer
layer number
Definition: SVDCoGTimeCalibrationsMonitorModule.h:70
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::SVDCoGTimeCalibrationsMonitorModule::beginRun
virtual void beginRun() override
initialize the TTrees and check validities of payloads
Definition: SVDCoGTimeCalibrationsMonitorModule.cc:36
Belle2::SVDCoGTimeCalibrationsMonitorModule::b_run
TBranch * b_run
run number
Definition: SVDCoGTimeCalibrationsMonitorModule.h:68
Belle2::SVDCoGTimeCalibrationsMonitorModule::b_triggerBin
TBranch * b_triggerBin
trigger bin
Definition: SVDCoGTimeCalibrationsMonitorModule.h:73
Belle2::SVDCoGTimeCalibrationsMonitorModule::m_rootFilePtr
TFile * m_rootFilePtr
pointer at root file used for storing histograms
Definition: SVDCoGTimeCalibrationsMonitorModule.h:64
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::Module
Base class for Modules.
Definition: Module.h:74
Belle2::SVDCoGTimeCalibrationsMonitorModule::b_sensor
TBranch * b_sensor
sensor number
Definition: SVDCoGTimeCalibrationsMonitorModule.h:71
Belle2::SVDCoGTimeCalibrationsMonitorModule::m_alpha
float m_alpha
alpha
Definition: SVDCoGTimeCalibrationsMonitorModule.h:85
Belle2::SVDCoGTimeCalibrationsMonitorModule::m_ladder
int m_ladder
ladder number
Definition: SVDCoGTimeCalibrationsMonitorModule.h:81
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::SVDCoGTimeCalibrationsMonitorModule::b_side
TBranch * b_side
sensor side
Definition: SVDCoGTimeCalibrationsMonitorModule.h:72
Belle2::SVDCoGTimeCalibrations::getCorrectedTime
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.
Definition: SVDCoGTimeCalibrations.h:66
Belle2::SVDCoGTimeCalibrationsMonitorModule::m_CoGTimeCal
SVDCoGTimeCalibrations m_CoGTimeCal
CoG time calibrations payload.
Definition: SVDCoGTimeCalibrationsMonitorModule.h:90
Belle2::VXD::SensorInfoBase::SVD
@ SVD
SVD Sensor.
Definition: SensorInfoBase.h:45
Belle2::SVDCoGTimeCalibrationsMonitorModule::m_tree
TTree * m_tree
pointer at tree containing the mean and RMS of calibration constants
Definition: SVDCoGTimeCalibrationsMonitorModule.h:65
Belle2::SVDCoGTimeCalibrationsMonitorModule::m_side
int m_side
sensor side
Definition: SVDCoGTimeCalibrationsMonitorModule.h:83
Belle2::SVDCoGTimeCalibrationsMonitorModule::m_run
int m_run
run number
Definition: SVDCoGTimeCalibrationsMonitorModule.h:79
Belle2::VXD::GeoCache
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:41
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::SVDCoGTimeCalibrationsMonitorModule::m_beta
float m_beta
beta
Definition: SVDCoGTimeCalibrationsMonitorModule.h:86
Belle2::SVDCoGTimeCalibrationsMonitorModule::b_alpha
TBranch * b_alpha
alpha
Definition: SVDCoGTimeCalibrationsMonitorModule.h:74
Belle2::SVDCoGTimeCalibrationsMonitorModule::endRun
virtual void endRun() override
print the payloads uniqueID and write tree and histograms to the rootfile
Definition: SVDCoGTimeCalibrationsMonitorModule.cc:113
Belle2::SVDCoGTimeCalibrations::isValid
bool isValid()
returns true if the m_aDBObtPtr is valid in the requested IoV
Definition: SVDCoGTimeCalibrations.h:115
Belle2::SVDCoGTimeCalibrationsMonitorModule::b_beta
TBranch * b_beta
beta
Definition: SVDCoGTimeCalibrationsMonitorModule.h:75
Belle2::SVDCoGTimeCalibrationsMonitorModule::m_sensor
int m_sensor
sensor number
Definition: SVDCoGTimeCalibrationsMonitorModule.h:82
Belle2::SVDCoGTimeCalibrationsMonitorModule::m_triggerBin
int m_triggerBin
trigger bin
Definition: SVDCoGTimeCalibrationsMonitorModule.h:84
Belle2::SVDCoGTimeCalibrationsMonitorModule::m_layer
int m_layer
layer number
Definition: SVDCoGTimeCalibrationsMonitorModule.h:80