Belle II Software  release-06-02-00
PXDTrackClusterDQMModule.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 
9 #include <pxd/modules/pxdDQM/PXDTrackClusterDQMModule.h>
10 #include <pxd/unpacking/PXDMappingLookup.h>
11 #include <pxd/geometry/SensorInfo.h>
12 #include <TDirectory.h>
13 
14 using namespace Belle2;
15 using namespace Belle2::PXD;
16 
17 //-----------------------------------------------------------------
18 // Register the Module
19 //-----------------------------------------------------------------
20 REG_MODULE(PXDTrackClusterDQM)
21 
22 
23 //-----------------------------------------------------------------
24 // Implementation
25 //-----------------------------------------------------------------
26 
27 PXDTrackClusterDQMModule::PXDTrackClusterDQMModule() : HistoModule(), m_vxdGeometry(VXD::GeoCache::getInstance())
28 {
29  //Set module properties
30  setDescription("DQM for PXD Cluster matched to a Track");
31  setPropertyFlags(c_ParallelProcessingCertified);
32 
33  addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms will be placed",
34  std::string("PXDER"));
35  addParam("moreHistos", m_moreHistos, "Fill additional histograms (not for ereco)", false);
36  addParam("ASICHistos", m_ASICHistos, "Fill additional histograms ASIC combination", true);
37 }
38 
39 
40 //------------------------------------------------------------------
41 // Function to define histograms
42 //-----------------------------------------------------------------
43 
45 {
46  m_recoTracks.isOptional(m_RecoTracksStoreArrayName);
47  m_tracks.isOptional(m_TracksStoreArrayName);
48 
49  // Register histograms (calls back defineHisto)
50  REG_HISTOGRAM
51 
52 }
53 
55 {
56  // Create a separate histogram directories and cd into it.
57  TDirectory* oldDir = gDirectory;
58  if (m_histogramDirectoryName != "") {
59  oldDir->mkdir(m_histogramDirectoryName.c_str());// do not use return value with ->cd(), its ZERO if dir already exists
60  oldDir->cd(m_histogramDirectoryName.c_str());
61  }
62 
63  std::vector<VxdID> sensors = m_vxdGeometry.getListOfSensors();
64  for (VxdID& avxdid : sensors) {
65  VXD::SensorInfoBase info = m_vxdGeometry.getSensorInfo(avxdid);
66  if (info.getType() != VXD::SensorInfoBase::PXD) continue;
67  //Only interested in PXD sensors
68 
69  TString buff = (std::string)avxdid;
70  buff.ReplaceAll(".", "_");
71 
72  m_trackClusterCharge[avxdid] = new TH1F("PXD_Track_Cluster_Charge_" + buff, "PXD Track Cluster Charge " + buff + ";Charge/ADU;",
73  100, 0, 100);
74  if (m_moreHistos) {
75  m_trackClusterChargeUC[avxdid] = new TH1F("PXD_Track_Cluster_Charge_UC_" + buff,
76  "PXD Track Cluster Charge (uncorrected)" + buff + ";Charge/ADU;", 100, 0, 100);
77  }
78  if (m_ASICHistos) {
79  // for now, we only want to have this module in
80  if (avxdid == VxdID("1.5.1")) {
81  for (int s = 0; s < 6; s++) {
82  for (int d = 0; d < 4; d++) {
83  m_trackASICClusterCharge[avxdid][s][d] = new TH1F("PXD_Track_Cluster_Charge_" + buff + Form("_sw%d_dcd%d", s + 1, d + 1),
84  "PXD Track Cluster Charge " + buff + Form(" sw%d dcd%d ", s + 1, d + 1) + ";Charge/ADU;", 100, 0, 100);
85  }
86  }
87  }
88  }
89  }
90 
91  m_trackedClusters = new TH1F("PXD_Tracked_Clusters", "PXD_Tracked_Clusters", 64, 0, 64);
92 
93  for (auto i = 0; i < 64; i++) {
94  auto layer = (((i >> 5) & 0x1) + 1);
95  auto ladder = ((i >> 1) & 0xF);
96  auto sensor = ((i & 0x1) + 1);
97 
98  auto id = Belle2::VxdID(layer, ladder, sensor);
99  // Check if sensor exist
100  if (Belle2::VXD::GeoCache::getInstance().validSensorID(id)) {
101  m_vxd_to_dhe[id] = i;
102  }
103  }
104 
105  oldDir->cd();
106 
107 }
108 
110 {
111  for (auto& it : m_trackClusterCharge) if (it.second) it.second->Reset();
112  for (auto& it : m_trackClusterChargeUC) if (it.second) it.second->Reset();
113  if (m_trackedClusters) m_trackedClusters->Reset();
114  for (const auto& it1 : m_trackASICClusterCharge) {
115  for (const auto& it2 : it1.second) {
116  for (const auto& it3 : it2) {
117  if (it3) it3->Reset();
118  }
119  }
120  }
121 }
122 
123 
125 {
126  if (m_trackedClusters) m_trackedClusters->Fill(-1); // Underflow as event counter
127  for (const Track& track : m_tracks) {
131  RelationVector<RecoTrack> recoTrack = track.getRelationsTo<RecoTrack>(m_RecoTracksStoreArrayName);
132  if (!recoTrack.size()) continue;
133  RelationVector<PXDCluster> pxdClustersTrack = DataStore::getRelationsWithObj<PXDCluster>(recoTrack[0]);
134 
135  const TrackFitResult* tfr = track.getTrackFitResultWithClosestMass(Const::pion);
136  double correction = 1.0;
137  if (tfr) correction = sin(tfr->getMomentum().Theta());
138  for (auto& cluster : pxdClustersTrack) {
139  m_trackedClusters->Fill(m_vxd_to_dhe[cluster.getSensorID()]);
140  if (m_trackClusterChargeUC[cluster.getSensorID()]) m_trackClusterChargeUC[cluster.getSensorID()]->Fill(cluster.getCharge());
141  if (tfr && m_trackClusterCharge[cluster.getSensorID()]) m_trackClusterCharge[cluster.getSensorID()]->Fill(
142  cluster.getCharge()*correction);
143  if (m_ASICHistos && tfr) {
144  // for now, we only want to have this module in
145  if (cluster.getSensorID() == VxdID("1.5.1")) {
146  auto SensorInfo = dynamic_cast<const PXD::SensorInfo&>(VXD::GeoCache::get(cluster.getSensorID()));
147  auto d = PXDMappingLookup::getDCDID(SensorInfo.getUCellID(cluster.getU()), SensorInfo.getVCellID(cluster.getV()),
148  cluster.getSensorID());
149  auto s = PXDMappingLookup::getSWBID(SensorInfo.getVCellID(cluster.getV()));
150 
151  TH1F* h = nullptr;
152  try {
153  h = m_trackASICClusterCharge[cluster.getSensorID()].at(s - 1).at(d - 1);
154  } catch (...) {
155  }
156  if (h) h->Fill(cluster.getCharge()*correction);
157  }
158  }
159  }
160  }
161 }
static const ChargedStable pion
charged pion particle
Definition: Const.h:542
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
DQM of cluster in matched tracks.
void initialize() override final
Initialize.
void defineHisto() override final
Define histograms.
void beginRun() override final
Begin run.
static int getDCDID(const int u, const int v, const VxdID sensorID)
get ID of DCD for giving pixel, range: 1..4.
static int getSWBID(const int v)
get ID of SWB for giving pixel, range: 1..6.
Specific implementation of SensorInfo for PXD Sensors which provides additional pixel specific inform...
Definition: SensorInfo.h:23
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:76
Class for type safe access to objects that are referred to in relations.
size_t size() const
Get number of relations.
Values of the result of a track fit with a given particle hypothesis.
TVector3 getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
Class that bundles various TrackFitResults.
Definition: Track.h:25
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:213
static const SensorInfoBase & get(Belle2::VxdID id)
Return a reference to the SensorInfo of a given SensorID.
Definition: GeoCache.h:139
Base class to provide Sensor Information for PXD and SVD.
int getVCellID(double v, bool clamp=false) const
Return the corresponding pixel/strip ID of a given v coordinate.
int getUCellID(double u, double v=0, bool clamp=false) const
Return the corresponding pixel/strip ID of a given u coordinate.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Namespace to encapsulate code needed for simulation and reconstrucion of the PXD.
Abstract base class for different kinds of events.