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