Belle II Software prerelease-11-00-00a
SVDClusterAbsoluteTimeShifterCollectorModule.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/SVDClusterAbsoluteTimeShifterCollectorModule.h>
9
10using namespace Belle2;
11
12//-----------------------------------------------------------------
13// Register the Module
14//-----------------------------------------------------------------
15REG_MODULE(SVDClusterAbsoluteTimeShifterCollector);
16
17//-----------------------------------------------------------------
18// Implementation
19//-----------------------------------------------------------------
20
22{
23 //Set module properties
24
25 setDescription("Collector module used to create the histograms needed for the SVD 6-Sample CoG, 3-Sample CoG and 3-Sample ELS Time calibration. "
26 "This is for the absolute shift of the calibrated time.");
28
29 addParam("SVDClustersOnTrackPrefix", m_svdClustersOnTrackPrefix, "set prefix of the list of clusters on track",
31 addParam("TimeAlgorithms", m_timeAlgorithms, "set list of time algorithms", m_timeAlgorithms);
32 addParam("EventT0Name", m_eventT0Name, "Name of the EventT0 list", m_eventT0Name);
33}
34
36{
37
38 for (auto alg : m_timeAlgorithms) {
40 m_svdClustersOnTrack[alg].isRequired((m_svdClustersOnTrackPrefix + '_' + alg).data());
41 }
42 m_eventT0.isRequired(m_eventT0Name);
43
44 //temporary histogram for checks: adding cdceventt0 to the histos
45 TH1F* __CDCEventT0__ = new TH1F("hCDCEventT0_",
46 "CDC Event T0",
47 300, -150., 150.);
48 __CDCEventT0__->GetXaxis()->SetTitle("CDC Event T0 (ns)");
49 registerObject<TH1F>(__CDCEventT0__->GetName(), __CDCEventT0__);
50
51 // getting all the svd sensors
53 std::vector<Belle2::VxdID> allSensors;
54 for (auto layer : geoCache.getLayers(VXD::SensorInfoBase::SVD))
55 for (auto ladder : geoCache.getLadders(layer)) {
56 for (Belle2::VxdID sensor : geoCache.getSensors(ladder))
57 allSensors.push_back(sensor);
58 break;
59 }
60
61 for (auto alg : m_timeAlgorithms) {
62
63 // 12 bins for L1-6 U-V sides
64
65 TH2F* __hClsOnTrack__ = new TH2F(("hClsTimeOnTracks_" + alg).data(),
66 ("Cluster time for " + alg).data(),
67 300, -150., 150., 12, 1., 13.);
68 __hClsOnTrack__->GetXaxis()->SetTitle("Cluster Time (ns)");
69 __hClsOnTrack__->GetYaxis()->SetTitle("LayerSideID");
70 registerObject<TH2F>(__hClsOnTrack__->GetName(), __hClsOnTrack__);
71 }
72
73}
74
76{
77 for (auto alg : m_timeAlgorithms) {
78 getObjectPtr<TH2F>(("hClsTimeOnTracks_" + alg).data())->Reset();
79 }
80 getObjectPtr<TH1F>("hCDCEventT0_")->Reset();
81}
82
84{
85 float eventT0 = 0;
86 // Set the CDC event t0 value if it exists
87 if (m_eventT0->hasTemporaryEventT0(Const::EDetector::CDC)) {
88 const auto evtT0CDC = m_eventT0->getBestCDCTemporaryEventT0();
89 eventT0 = evtT0CDC->eventT0;
90 } else {return;}
91
92 // Fill histograms clusters on tracks
93 for (auto alg : m_timeAlgorithms) {
94
95 if (!m_svdClustersOnTrack[alg].isValid()) {
96 B2WARNING("!!!! List is not Valid: isValid() = " << m_svdClustersOnTrack[alg].isValid());
97 return;
98 }
99 for (const auto& svdCluster : m_svdClustersOnTrack[alg]) {
100
101 // get cluster time
102 float clusterTime = svdCluster.getClsTime();
103
104 int side = svdCluster.isUCluster() ? 0 : 1;
105 int layerNumber = svdCluster.getSensorID().getLayerNumber();
106 int LayerSideID = 2 * layerNumber - side ; // 1 to 12 something similar to sensorBin
107 double LayerSideBinCenter = getObjectPtr<TH2F>("hClsTimeOnTracks_" + alg)->GetYaxis()->GetBinCenter(LayerSideID);
108
109 getObjectPtr<TH2F>(("hClsTimeOnTracks_" + alg).data())->Fill(clusterTime - eventT0, LayerSideBinCenter);
110
111 }
112 } // loop over alg
113
114 getObjectPtr<TH1F>("hCDCEventT0_")->Fill(eventT0);
115}
void registerObject(const std::string &name, T *obj)
Register object with a name, takes ownership, do not access the pointer beyond prepare()
T * getObjectPtr(const std::string &name)
Calls the CalibObjManager to get the requested stored collector data.
CalibrationCollectorModule()
Constructor. Sets the default prefix for calibration dataobjects.
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
std::string m_svdClustersOnTrackPrefix
Name Prefix of the SVDClustersOnTracks store array.
std::string m_eventT0Name
Name of the EventT0 store object pointer used as parameter of the module.
std::vector< std::string > m_timeAlgorithms
List of time algorithms to calibrate.
std::map< TString, StoreArray< SVDCluster > > m_svdClustersOnTrack
SVDClusters store array.
Accessor to arrays stored in the data store.
Definition StoreArray.h:113
Class to facilitate easy access to sensor information of the VXD like coordinate transformations or p...
Definition GeoCache.h:38
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:32
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:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition Module.h:649
Abstract base class for different kinds of events.