Belle II Software  release-08-01-10
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 
22 SVDTimeCalibrationCollectorModule::SVDTimeCalibrationCollectorModule() : CalibrationCollectorModule()
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");
28 
29  addParam("SVDClustersFromTracksName", m_svdClustersOnTracks, "Name of the SVDClusters list", m_svdClustersOnTracks);
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  addParam("RawTimeIoVMin", m_minRawTimeForIoV,
34  "Minimum value of the raw time distribution used to determine whether change IoV or not", m_minRawTimeForIoV);
35  addParam("RawTimeIoVMax", m_maxRawTimeForIoV,
36  "Maximum value of the raw time distribution used to determine whether change IoV or not", m_maxRawTimeForIoV);
37 }
38 
40 {
41 
42  m_hEventT0FromCDC = new TH1F("hEventT0FromCDC", "EventT0FromCDC", 200, -100, 100);
43  registerObject<TH1F>("hEventT0FromCDC", m_hEventT0FromCDC);
44  m_hEventT0FromCDCSync = new TH1F("hEventT0FromCDCSync", "EventT0FromCDCSync", 200, -100, 100);
45  registerObject<TH1F>("hEventT0FromCDCSync", m_hEventT0FromCDCSync);
46  m_hRawTimeL3V = new TH1F("hRawTimeL3V", "RawTimeL3V", 150, 0, 150);
47  registerObject<TH1F>("hRawTimeL3V", m_hRawTimeL3V);
48  m_hRawTimeL3VFullRange = new TH1F("hRawTimeL3VFullRange", "RawTimeL3V full range", 400, -150, 250);
49  registerObject<TH1F>("hRawTimeL3VFullRange", m_hRawTimeL3VFullRange);
50 
52  m_eventT0.isRequired(m_eventTime);
53  // m_svdEI.isRequired(m_svdEventInfo);
54 
56  std::vector<Belle2::VxdID> allSensors;
57  for (auto layer : geoCache.getLayers(VXD::SensorInfoBase::SVD))
58  for (auto ladder : geoCache.getLadders(layer))
59  for (Belle2::VxdID sensor : geoCache.getSensors(ladder))
60  allSensors.push_back(sensor);
61 
62  int numberOfSensorBin = 2 * int(allSensors.size());
63  B2INFO("Number of SensorBin: " << numberOfSensorBin);
64 
65  TH3F* __hEventT0vsCoG__ = new TH3F("__hEventT0vsCoG__", "EventT0Sync vs rawTime",
66  int(200 / m_rawCoGBinWidth), -100, 100, 60, -100, 20,
67  numberOfSensorBin, + 0.5, numberOfSensorBin + 0.5);
68  TH2F* __hEventT0__ = new TH2F("__hEventT0__", "EventT0Sync",
69  100, -100, 100,
70  numberOfSensorBin, + 0.5, numberOfSensorBin + 0.5);
71  TH2F* __hEventT0NoSync__ = new TH2F("__hEventT0NoSync__", "EventT0NoSync",
72  100, -100, 100,
73  numberOfSensorBin, + 0.5, numberOfSensorBin + 0.5);
74  TH1F* __hBinToSensorMap__ = new TH1F("__hBinToSensorMap__", "__BinToSensorMap__",
75  numberOfSensorBin, + 0.5, numberOfSensorBin + 0.5);
76  __hEventT0vsCoG__->GetZaxis()->SetTitle("Sensor");
77  __hEventT0vsCoG__->GetYaxis()->SetTitle("EventT0Sync (ns)");
78  __hEventT0vsCoG__->GetXaxis()->SetTitle("raw_time (ns)");
79  __hEventT0__->GetYaxis()->SetTitle("Sensor");
80  __hEventT0__->GetXaxis()->SetTitle("event_t0 (ns)");
81  __hEventT0NoSync__->GetYaxis()->SetTitle("sensor");
82  __hEventT0NoSync__->GetXaxis()->SetTitle("event_t0 (ns)");
83 
84  int tmpBinCnt = 0;
85  for (auto sensor : allSensors) {
86  for (auto view : {'U', 'V'}) {
87  tmpBinCnt++;
88  TString binLabel = TString::Format("L%iL%iS%i%c",
89  sensor.getLayerNumber(),
90  sensor.getLadderNumber(),
91  sensor.getSensorNumber(),
92  view);
93  __hBinToSensorMap__->GetXaxis()->SetBinLabel(tmpBinCnt, binLabel);
94  }
95  }
96  registerObject<TH3F>(__hEventT0vsCoG__->GetName(), __hEventT0vsCoG__);
97  registerObject<TH2F>(__hEventT0__->GetName(), __hEventT0__);
98  registerObject<TH2F>(__hEventT0NoSync__->GetName(), __hEventT0NoSync__);
99  registerObject<TH1F>(__hBinToSensorMap__->GetName(), __hBinToSensorMap__);
100 }
101 
103 {
104  getObjectPtr<TH1F>("hEventT0FromCDC")->Reset();
105  getObjectPtr<TH1F>("hEventT0FromCDCSync")->Reset();
106  getObjectPtr<TH1F>("hRawTimeL3V")->Reset();
107  getObjectPtr<TH1F>("hRawTimeL3VFullRange")->Reset();
108  getObjectPtr<TH3F>("__hEventT0vsCoG__")->Reset();
109  getObjectPtr<TH2F>("__hEventT0__")->Reset();
110  getObjectPtr<TH2F>("__hEventT0NoSync__")->Reset();
111  getObjectPtr<TH1F>("__hBinToSensorMap__")->Reset();
112 }
113 
115 {
116  float eventT0 = 0;
117  // Set the CDC event t0 value if it exists
118  if (m_eventT0->hasTemporaryEventT0(Const::EDetector::CDC)) {
119  auto evtT0CDC = m_eventT0->getBestCDCTemporaryEventT0();
120  // Set the CDC event t0 value for filling into the histogram
121  // The most accurate CDC event t0 value is the last one in the list.
122  eventT0 = evtT0CDC->eventT0;
123  getObjectPtr<TH1F>("hEventT0FromCDC")->Fill(eventT0);
124  } else {return;}
125 
126  if (!m_svdClsOnTrk.isValid()) {
127  B2WARNING("!!!! File is not Valid: isValid() = " << m_svdClsOnTrk.isValid());
128  return;
129  }
130 
131  //get SVDEventInfo
132  StoreObjPtr<SVDEventInfo> temp_eventinfo("SVDEventInfo");
133  std::string m_svdEventInfoName = "SVDEventInfo";
134  if (!temp_eventinfo.isOptional("SVDEventInfo"))
135  m_svdEventInfoName = "SVDEventInfoSim";
136  StoreObjPtr<SVDEventInfo> eventinfo(m_svdEventInfoName);
137  if (!eventinfo) B2ERROR("No SVDEventInfo!");
138  eventinfo->setAPVClock(m_hwClock);
139 
140  for (const auto& svdCluster : m_svdClsOnTrk) {
141  // get cluster time
142  float clTime_ftsw = svdCluster.getClsTime();
143 
144  //remove firstFrame and triggerBin correction applied in the clusterizer AND relative shift 3/6 mixed sample DAQ
145  float clTime = eventinfo->getTimeInSVDReference(clTime_ftsw, svdCluster.getFirstFrame());
146 
147  //get cluster side
148  int side = svdCluster.isUCluster();
149 
150  //get layer
151  short unsigned int layer = svdCluster.getSensorID().getLayerNumber();
152 
153  //fill histograms only if EventT0 is there
154  if (m_eventT0->hasTemporaryEventT0(Const::EDetector::CDC)) {
155 
156  float eventT0Sync = eventinfo->getTimeInSVDReference(eventT0, svdCluster.getFirstFrame());
157 
158  TString binLabel = TString::Format("L%iL%iS%i%c",
159  svdCluster.getSensorID().getLayerNumber(),
160  svdCluster.getSensorID().getLadderNumber(),
161  svdCluster.getSensorID().getSensorNumber(),
162  side ? 'U' : 'V');
163  int sensorBin = getObjectPtr<TH1F>("__hBinToSensorMap__")->GetXaxis()->FindBin(binLabel.Data());
164  double sensorBinCenter = getObjectPtr<TH1F>("__hBinToSensorMap__")->GetXaxis()->GetBinCenter(sensorBin);
165  getObjectPtr<TH3F>("__hEventT0vsCoG__")->Fill(clTime, eventT0Sync, sensorBinCenter);
166  getObjectPtr<TH2F>("__hEventT0__")->Fill(eventT0Sync, sensorBinCenter);
167  getObjectPtr<TH2F>("__hEventT0NoSync__")->Fill(eventT0, sensorBinCenter);
168 
169  getObjectPtr<TH1F>("hEventT0FromCDCSync")->Fill(eventT0Sync);
170  if (layer == 3 && side == 0) {
171  if (clTime_ftsw >= m_minRawTimeForIoV && clTime_ftsw <= m_maxRawTimeForIoV) {
172  getObjectPtr<TH1F>("hRawTimeL3V")->Fill(clTime_ftsw);
173  }
174  getObjectPtr<TH1F>("hRawTimeL3VFullRange")->Fill(clTime_ftsw);
175  }
176  }
177  };
178 }
Calibration collector module base class.
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
StoreArray< SVDCluster > m_svdClsOnTrk
SVDClusters store array.
double m_rawCoGBinWidth
Raw_CoG Bin Width [ns] for 2D-histogram.
void prepare() override final
Initialize the module.
StoreObjPtr< EventT0 > m_eventT0
EventT0 store object pointer.
TH1F * m_hEventT0FromCDC
Distribution of EventT0 reconstructed by the CDC for all sensos/side.
double m_maxRawTimeForIoV
Maxmum value of the raw time distribution used to determine whether change IoV or not.
double m_minRawTimeForIoV
Minimum value of the raw time distribution used to determine whether change IoV or not.
TH1F * m_hEventT0FromCDCSync
Distribution of EventT0 reconstructed by the CDC and synchronized for all sensos/side.
TH1F * m_hRawTimeL3V
Raw time distribution of layer3 V-side for IoV determination.
DBObjPtr< HardwareClockSettings > m_hwClock
systems clock
TH1F * m_hRawTimeL3VFullRange
Raw time distribution of layer3 V-side.
void startRun() override final
Called when entering a new run.
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
bool isOptional(const std::string &name="")
Tell the DataStore about an optional input.
bool isValid() const
Check wether the array was registered.
Definition: StoreArray.h:288
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
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
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.