Belle II Software  release-06-02-00
SVDTimeCalibrationCollectorModule.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/svdTimeCalibrationCollector/SVDTimeCalibrationCollectorModule.h>
9 
10 using namespace std;
11 using namespace Belle2;
12 
13 //-----------------------------------------------------------------
14 // Register the Module
15 //-----------------------------------------------------------------
16 REG_MODULE(SVDTimeCalibrationCollector)
17 
18 //-----------------------------------------------------------------
19 // Implementation
20 //-----------------------------------------------------------------
21 
23 {
24  //Set module properties
25 
26  setDescription("Collector module used to create the histograms needed for the SVD 6-Sample CoG, 3-Sample CoG and 3-Sample ELS Time calibration");
27  setPropertyFlags(c_ParallelProcessingCertified);
28 
29  addParam("SVDClustersFromTracksName", m_svdClusters, "Name of the SVDClusters list", m_svdClusters);
30  addParam("EventT0Name", m_eventTime, "Name of the EventT0 list", m_eventTime);
31  addParam("SVDEventInfoName", m_svdEventInfo, "Name of the SVDEventInfo list", m_svdEventInfo);
32  addParam("RawCoGBinWidth", m_rawCoGBinWidth, "Bin Width [ns] for raw CoG time", m_rawCoGBinWidth);
33 }
34 
35 void SVDTimeCalibrationCollectorModule::prepare()
36 {
37  TH2F hEventT0vsCoG("eventT0vsCoG__L@layerL@ladderS@sensor@view",
38  "EventT0Sync vs rawTime in @layer.@ladder.@sensor @view/@side",
39  int(200 / m_rawCoGBinWidth), -100, 100, 100, -100, 100);
40  hEventT0vsCoG.GetYaxis()->SetTitle("EventT0Sync (ns)");
41  hEventT0vsCoG.GetXaxis()->SetTitle("raw_time (ns)");
42  m_hEventT0vsCoG = new SVDHistograms<TH2F>(hEventT0vsCoG);
43 
44  TH1F hEventT0("eventT0__L@layerL@ladderS@sensor@view",
45  "EventT0Sync in @layer.@ladder.@sensor @view/@side",
46  100, -100, 100);
47  hEventT0.GetXaxis()->SetTitle("event_t0 (ns)");
48  m_hEventT0 = new SVDHistograms<TH1F>(hEventT0);
49 
50  TH1F hEventT0NoSync("eventT0nosync__L@layerL@ladderS@sensor@view",
51  "EventT0NoSync in @layer.@ladder.@sensor @view/@side",
52  100, -100, 100);
53  hEventT0NoSync.GetXaxis()->SetTitle("event_t0 (ns)");
54  m_hEventT0nosync = new SVDHistograms<TH1F>(hEventT0NoSync);
55 
56  m_hEventT0FromCDC = new TH1F("hEventT0FromCDC", "EventT0FromCDC", 200, -100, 100);
57  registerObject<TH1F>("hEventT0FromCDC", m_hEventT0FromCDC);
58  m_hEventT0FromCDCSync = new TH1F("hEventT0FromCDCSync", "EventT0FromCDCSync", 200, -100, 100);
59  registerObject<TH1F>("hEventT0FromCDCSync", m_hEventT0FromCDCSync);
60  m_hRawTimeL3V = new TH1F("hRawTimeL3V", "RawCoGTimeL3V", 300, -150, 150);
61  registerObject<TH1F>("hRawTimeL3V", m_hRawTimeL3V);
62 
63  m_svdCls.isRequired(m_svdClusters);
64  m_eventT0.isRequired(m_eventTime);
65  m_svdEI.isRequired(m_svdEventInfo);
66 
67  VXD::GeoCache& geoCache = VXD::GeoCache::getInstance();
68 
69  for (auto layer : geoCache.getLayers(VXD::SensorInfoBase::SVD)) {
70  for (auto ladder : geoCache.getLadders(layer)) {
71  for (Belle2::VxdID sensor : geoCache.getSensors(ladder)) {
72  for (int view = SVDHistograms<TH2F>::VIndex ; view < SVDHistograms<TH2F>::UIndex + 1; view++) {
73  registerObject<TH2F>(m_hEventT0vsCoG->getHistogram(sensor, view)->GetName(), m_hEventT0vsCoG->getHistogram(sensor, view));
74  registerObject<TH1F>(m_hEventT0->getHistogram(sensor, view)->GetName(), m_hEventT0->getHistogram(sensor, view));
75  registerObject<TH1F>(m_hEventT0nosync->getHistogram(sensor, view)->GetName(), m_hEventT0nosync->getHistogram(sensor, view));
76  }
77  }
78  }
79  }
80 }
81 
82 void SVDTimeCalibrationCollectorModule::startRun()
83 {
84 
85  VXD::GeoCache& geoCache = VXD::GeoCache::getInstance();
86 
87  for (auto layer : geoCache.getLayers(VXD::SensorInfoBase::SVD)) {
88  for (auto ladder : geoCache.getLadders(layer)) {
89  for (Belle2::VxdID sensor : geoCache.getSensors(ladder)) {
90  for (int view = SVDHistograms<TH2F>::VIndex ; view < SVDHistograms<TH2F>::UIndex + 1; view++) {
91  // std::string s = std::string(sensor);
92  // std::string v = std::to_string(view);
93  // std::string name = string("eventT0vsCog_")+s+string("_")+v;
94  // registerObject<TH2F>(name.c_str(),m_hEventT0vsCoG->getHistogram(sensor, view));
95  getObjectPtr<TH2F>(m_hEventT0vsCoG->getHistogram(sensor, view)->GetName())->Reset();
96  getObjectPtr<TH1F>(m_hEventT0->getHistogram(sensor, view)->GetName())->Reset();
97  getObjectPtr<TH1F>(m_hEventT0nosync->getHistogram(sensor, view)->GetName())->Reset();
98  }
99  }
100  }
101  }
102 }
103 
104 void SVDTimeCalibrationCollectorModule::collect()
105 {
106  float eventT0 = 0;
107  // Set the CDC event t0 value if it exists
108  if (m_eventT0->hasTemporaryEventT0(Const::EDetector::CDC)) {
109  auto evtT0List_CDC = m_eventT0->getTemporaryEventT0s(Const::EDetector::CDC);
110  // Set the CDC event t0 value for filling into the histogram
111  // The most accurate CDC event t0 value is the last one in the list.
112  eventT0 = evtT0List_CDC.back().eventT0;
113  getObjectPtr<TH1F>("hEventT0FromCDC")->Fill(eventT0);
114  } else {return;}
115 
116  if (!m_svdCls.isValid()) {
117  B2WARNING("!!!! File is not Valid: isValid() = " << m_svdCls.isValid());
118  return;
119  }
120 
121  //get SVDEventInfo
122  StoreObjPtr<SVDEventInfo> temp_eventinfo("SVDEventInfo");
123  std::string m_svdEventInfoName = "SVDEventInfo";
124  if (!temp_eventinfo.isOptional("SVDEventInfo"))
125  m_svdEventInfoName = "SVDEventInfoSim";
126  StoreObjPtr<SVDEventInfo> eventinfo(m_svdEventInfoName);
127  if (!eventinfo) B2ERROR("No SVDEventInfo!");
128 
129  for (int cl = 0 ; cl < m_svdCls.getEntries(); cl++) {
130  // get cluster time
131  float clTime_ftsw = m_svdCls[cl]->getClsTime();
132 
133  //remove firstFrame and triggerBin correction applied in the clusterizer AND relative shift 3/6 mixed sample DAQ
134  float clTime = eventinfo->getTimeInSVDReference(clTime_ftsw, m_svdCls[cl]->getFirstFrame());
135 
136  //get cluster side
137  int side = m_svdCls[cl]->isUCluster();
138 
139  //get VxdID
140  VxdID::baseType theVxdID = (VxdID::baseType)m_svdCls[cl]->getSensorID();
141  short unsigned int layer = m_svdCls[cl]->getSensorID().getLayerNumber();
142 
143  //fill histograms only if EventT0 is there
144  if (m_eventT0->hasTemporaryEventT0(Const::EDetector::CDC)) {
145 
146  float eventT0Sync = eventinfo->getTimeInSVDReference(eventT0, m_svdCls[cl]->getFirstFrame());
147 
148  getObjectPtr<TH2F>(m_hEventT0vsCoG->getHistogram(theVxdID, side)->GetName())->Fill(clTime, eventT0Sync);
149  getObjectPtr<TH1F>(m_hEventT0->getHistogram(theVxdID, side)->GetName())->Fill(eventT0Sync);
150  getObjectPtr<TH1F>(m_hEventT0nosync->getHistogram(theVxdID, side)->GetName())->Fill(eventT0);
151  getObjectPtr<TH1F>("hEventT0FromCDCSync")->Fill(eventT0Sync);
152  if (layer == 3 && side == 0) {getObjectPtr<TH1F>("hRawTimeL3V")->Fill(clTime_ftsw);}
153  }
154  };
155 }
Calibration collector module base class.
Collector module used to create the histograms needed for the SVD CoG-Time calibration.
bool isOptional(const std::string &name="")
Tell the DataStore about an optional input.
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:95
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:39
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
unsigned short baseType
The base integer type for VxdID.
Definition: VxdID.h:36
#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.