Belle II Software  release-05-01-25
PXDDQMEfficiencyModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2013 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Bjoern Spruck, Thomas, Uwe *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <pxd/modules/pxdDQM/PXDDQMEfficiencyModule.h>
12 #include <tracking/dataobjects/ROIid.h>
13 
14 #include <pxd/reconstruction/PXDPixelMasker.h>
15 #include <mdst/dataobjects/Track.h>
16 #include <framework/gearbox/Const.h>
17 
18 #include "TDirectory.h"
19 #include "TMatrixDSym.h"
20 using namespace Belle2;
21 
22 //-----------------------------------------------------------------
23 // Register the Module
24 //-----------------------------------------------------------------
25 REG_MODULE(PXDDQMEfficiency)
26 
27 //-----------------------------------------------------------------
28 // Implementation
29 //-----------------------------------------------------------------
30 
31 PXDDQMEfficiencyModule::PXDDQMEfficiencyModule() : HistoModule(), m_vxdGeometry(VXD::GeoCache::getInstance())
32 {
33  // Set module properties
34  setDescription("Create basic histograms for PXD efficiency");
35 
36  // What exactly is needed for this to be true?
37  setPropertyFlags(c_ParallelProcessingCertified);
38 
39  // Parameter definitions
40  addParam("pxdClustersName", m_pxdClustersName, "name of StoreArray with PXD cluster", std::string(""));
41  addParam("tracksName", m_tracksName, "name of StoreArray with RecoTracks", std::string(""));
42  addParam("ROIsName", m_ROIsName, "name of the list of HLT ROIs, if available in output", std::string(""));
43  addParam("PXDInterceptListName", m_PXDInterceptListName, "name of the list of interceptions", std::string(""));
44  addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms will be placed",
45  std::string("PXDEFF"));
46  addParam("binsU", m_u_bins, "histogram bins in u direction", int(4));
47  addParam("binsV", m_v_bins, "histogram bins in v direction", int(6));
48  addParam("distCut", m_distcut, "max distance in [cm] for cluster to be counted to a track", double(0.0500));
49  addParam("pCut", m_pcut, "Set a cut on the track p-value", double(1e-20));
50  addParam("requireROIs", m_requireROIs, "require tracks to lie inside a ROI", bool(false));
51  addParam("useAlignment", m_useAlignment, "if true the alignment will be used", bool(true));
52  addParam("maskDeadPixels", m_maskDeadPixels, "Do not consider tracks going through known dead or hot pixels for the efficiency",
53  bool(false));
54  addParam("minSVDHits", m_minSVDHits, "Number of SVD hits required in a track to be considered", 0u);
55  addParam("momCut", m_momCut, "Set a cut on the track momentum in GeV/c, 0 disables", double(0));
56  addParam("pTCut", m_pTCut, "Set a cut on the track pT in GeV/c, 0 disables", double(1));
57  addParam("cutBorders", m_cutBorders, "Do not use tracks near the borders of the sensor", bool(true));
58  addParam("maskedDistance", m_maskedDistance, "Distance inside which no masked pixel or sensor border is allowed", int(10));
59  addParam("trackUFactorDistCut", m_uFactor, "Set a cut on u error of track (factor*err<dist), 0 disables", double(2.0));
60  addParam("trackVFactorDistCut", m_vFactor, "Set a cut on v error of track (factor*err<dist), 0 disables", double(2.0));
61  addParam("z0minCut", m_z0minCut, "Set a cut z0 minimum in cm (large negativ value eg -9999 disables)", double(-1));
62  addParam("z0maxCut", m_z0maxCut, "Set a cut z0 maximum in cm (large positiv value eg 9999 disables)", double(1));
63  addParam("d0Cut", m_d0Cut, "Set a cut abs(d0) in cm (and negativ value eg -9999 disables)", double(0.5));
64  addParam("verboseHistos", m_verboseHistos, "Add more verbose histograms for cuts (not for ereoc)", bool(false));
65 }
66 
67 
69 {
70  //calls the define histogram function
71  REG_HISTOGRAM;
72 
73  //register the required arrays
74  //Register as optional so validation for cases where they are not available still succeeds, but module will not do any meaningful work without them
75  m_pxdclusters.isOptional(m_pxdClustersName);
76  m_tracks.isOptional(m_tracksName);
77  m_ROIs.isOptional(m_ROIsName);
79 }
80 
82 {
83  for (auto& h : m_h_track_hits) if (h.second) h.second->Reset();
84  for (auto& h : m_h_matched_cluster) if (h.second) h.second->Reset();
85  for (auto& h : m_h_p) if (h.second) h.second->Reset();
86  for (auto& h : m_h_pt) if (h.second) h.second->Reset();
87  for (auto& h : m_h_su) if (h.second) h.second->Reset();
88  for (auto& h : m_h_sv) if (h.second) h.second->Reset();
89  for (auto& h : m_h_p2) if (h.second) h.second->Reset();
90  for (auto& h : m_h_pt2) if (h.second) h.second->Reset();
91  for (auto& h : m_h_su2) if (h.second) h.second->Reset();
92  for (auto& h : m_h_sv2) if (h.second) h.second->Reset();
93 }
94 
96 {
97  if (!m_pxdclusters.isValid()) {
98  B2INFO("PXDClusters array is missing, no efficiencies");
99  return;
100  }
101  if (!m_tracks.isValid()) {
102  B2INFO("RecoTrack array is missing, no efficiencies");
103  return;
104  }
105  if (!m_ROIs.isValid() && m_requireROIs) {
106  B2INFO("ROI array is missing but required hits in ROIs, aborting");
107  return;
108  }
109  if (!m_intercepts.isValid()) {
110  B2INFO("Intercept array is missing, no efficiencies");
111  return;
112  }
113 
114 
115  for (auto& a_track : m_tracks) {
116 
117  //If fit failed assume position pointed to is useless anyway
118  if (!a_track.wasFitSuccessful()) continue;
119 
120  if (a_track.getNumberOfSVDHits() < m_minSVDHits) continue;
121 
122  RelationVector<PXDIntercept> interceptList = a_track.getRelationsTo<PXDIntercept>(m_PXDInterceptListName);
123  if (!interceptList.size()) continue;
124 
125  const genfit::FitStatus* fitstatus = a_track.getTrackFitStatus();
126  if (fitstatus->getPVal() < m_pcut) continue;
127 
128  genfit::MeasuredStateOnPlane trackstate;
129  trackstate = a_track.getMeasuredStateOnPlaneFromFirstHit();
130  if (trackstate.getMom().Mag() < m_momCut) continue;
131  if (trackstate.getMom().Pt() < m_pTCut) continue;
132 
133  auto ptr = a_track.getRelated<Track>("Tracks");
134 
135  if (!ptr) {
136  B2ERROR("expect a track for fitted recotracks");
137  continue;
138  }
139  auto ptr2 = ptr->getTrackFitResultWithClosestMass(Const::pion);
140  if (!ptr2) {
141  B2ERROR("expect a track fit result for mass");
142  continue;
143  }
144 
145  // Vertex cut
146  if (ptr2->getZ0() < m_z0minCut || ptr2->getZ0() > m_z0maxCut || fabs(ptr2->getD0()) > m_d0Cut) continue;
147 
148  //loop over all PXD sensors to get the intersections
149  std::vector<VxdID> sensors = m_vxdGeometry.getListOfSensors();
150  for (auto intercept : interceptList) {
151  auto const aVxdID = intercept.getSensorID();
152  auto& info = m_vxdGeometry.getSensorInfo(aVxdID);
153  //Search for intersections of the track with all PXD layers
154  //Traditional (aka the person before did it like this) method
155  //If there is a way to find out sensors crossed by a track directly, that would most likely be faster
156 
157  //true = track intersects current sensor
158  double sigu(-9999);
159  double sigv(-9999);
160  {
161  if (m_verboseHistos) {
162  if (m_h_p[aVxdID]) m_h_p[aVxdID]->Fill(trackstate.getMom().Mag());
163  if (m_h_pt[aVxdID]) m_h_pt[aVxdID]->Fill(trackstate.getMom().Pt());
164  if (m_h_su[aVxdID]) m_h_su[aVxdID]->Fill(sigu);
165  if (m_h_sv[aVxdID]) m_h_sv[aVxdID]->Fill(sigv);
166  }
167  if (m_uFactor * sigu > m_distcut) continue; // Error ufak*SigmaU > cut
168  if (m_vFactor * sigv > m_distcut) continue; // Error vfak*SigmaV > cut
169 
170  double u_fit = intercept.getCoorU();
171  double v_fit = intercept.getCoorV();
172 
173  int ucell_fit = info.getUCellID(u_fit); // check wie overflow!!!
174  int vcell_fit = info.getVCellID(v_fit); // Check wie overflow
175 
176  if (m_cutBorders && isCloseToBorder(ucell_fit, vcell_fit, m_maskedDistance)) {
177  continue;
178  }
179 
180  if (m_maskDeadPixels && isDeadPixelClose(ucell_fit, vcell_fit, m_maskedDistance, aVxdID)) {
181  continue;
182  }
183 
184  if (m_requireROIs) {
185  //Check if the intersection is inside a ROI
186  //If not, even if measured the cluster was thrown away->Not PXD's fault
187  bool fitInsideROI = false;
188  for (auto& roit : m_ROIs) {
189  if (aVxdID != roit.getSensorID()) {
190  continue; //ROI on other sensor
191  }
192 
193  if (ucell_fit < roit.getMaxUid()
194  && ucell_fit > roit.getMinUid()
195  && vcell_fit < roit.getMaxVid()
196  && vcell_fit > roit.getMinVid()) {
197  fitInsideROI = true;
198  }
199  }
200  if (!fitInsideROI) {
201  continue;//Hit wouldn't have been recorded
202  }
203  }
204 
205  //This track should be on the sensor
206  m_h_track_hits[aVxdID]->Fill(ucell_fit, vcell_fit);
207 
208  //Now check if the sensor measured a hit here
209 
210  int bestcluster = findClosestCluster(aVxdID, TVector3(u_fit, v_fit, 0));
211  if (bestcluster >= 0) {
212  double u_clus = m_pxdclusters[bestcluster]->getU();
213  double v_clus = m_pxdclusters[bestcluster]->getV();
214 
215  //is the closest cluster close enough to the track to count as measured?
216  TVector3 dist_clus(u_fit - u_clus, v_fit - v_clus, 0);
217  if (dist_clus.Mag() <= m_distcut) {
218  m_h_matched_cluster[aVxdID]->Fill(ucell_fit, vcell_fit);
219  if (m_verboseHistos) {
220  if (m_h_p2[aVxdID]) m_h_p2[aVxdID]->Fill(trackstate.getMom().Mag());
221  if (m_h_pt2[aVxdID]) m_h_pt2[aVxdID]->Fill(trackstate.getMom().Pt());
222  if (m_h_su2[aVxdID]) m_h_su2[aVxdID]->Fill(sigu);
223  if (m_h_sv2[aVxdID]) m_h_sv2[aVxdID]->Fill(sigv);
224  }
225  }
226  }
227  }
228  }
229  }
230 }
231 
233 {
234  // Create a separate histogram directory and cd into it.
235  TDirectory* oldDir = gDirectory;
236  if (m_histogramDirectoryName != "") {
237  oldDir->mkdir(m_histogramDirectoryName.c_str());
238  oldDir->cd(m_histogramDirectoryName.c_str());
239  }
240 
241  std::vector<VxdID> sensors = m_vxdGeometry.getListOfSensors();
242  for (VxdID& avxdid : sensors) {
244  if (info.getType() != VXD::SensorInfoBase::PXD) continue;
245  //Only interested in PXD sensors
246 
247  TString buff = (std::string)avxdid;
248  buff.ReplaceAll(".", "_");
249 
250  int nu = info.getUCells();
251  int nv = info.getVCells();
252 
253  //nu + 1,nv + 1 Bin numbers when using one bin per pixel
254  m_h_track_hits[avxdid] = new TH2F("track_hits_" + buff, "tracks through sensor " + buff,
255  m_u_bins, -0.5, nu - 0.5, m_v_bins, -0.5, nv - 0.5);
256  m_h_matched_cluster[avxdid] = new TH2F("matched_cluster_" + buff, "clusters matched to track intersections " + buff,
257  m_u_bins, -0.5, nu - 0.5, m_v_bins, -0.5, nv - 0.5);
258 
259  if (m_verboseHistos) {
260  m_h_p[avxdid] = new TH1F("p_" + buff, "p " + buff, 100, 0, 10);
261  m_h_pt[avxdid] = new TH1F("pt_" + buff, "pt " + buff, 100, 0, 10);
262  m_h_su[avxdid] = new TH1F("su_" + buff, "su " + buff, 1000, 0, 1);
263  m_h_sv[avxdid] = new TH1F("sv_" + buff, "sv " + buff, 1000, 0, 1);
264  m_h_p2[avxdid] = new TH1F("p2_" + buff, "p2 " + buff, 100, 0, 10);
265  m_h_pt2[avxdid] = new TH1F("pt2_" + buff, "pt2 " + buff, 100, 0, 10);
266  m_h_su2[avxdid] = new TH1F("su2_" + buff, "su2 " + buff, 1000, 0, 1);
267  m_h_sv2[avxdid] = new TH1F("sv2_" + buff, "sv2 " + buff, 1000, 0, 1);
268  }
269  }
270  // cd back to root directory
271  oldDir->cd();
272 }
273 
274 
275 int
276 PXDDQMEfficiencyModule::findClosestCluster(const VxdID& avxdid, TVector3 intersection)
277 {
278  int closest = -1;
279  double mindist = 999999999999; //definitely outside of the sensor
280 
282 
283  //loop the clusters
284  for (int iclus = 0; iclus < m_pxdclusters.getEntries(); iclus++) {
285  //Do not consider as different if only segment differs!
286  //As of this writing segment is never filled for clusters, but just to be sure
287  VxdID clusterID = m_pxdclusters[iclus]->getSensorID();
288  if (avxdid.getLayerNumber() != clusterID.getLayerNumber() ||
289  avxdid.getLadderNumber() != clusterID.getLadderNumber() ||
290  avxdid.getSensorNumber() != clusterID.getSensorNumber()) {
291  continue;
292  }
293  //only cluster on the correct sensor and direction should survive
294 
295  double u = m_pxdclusters[iclus]->getU();
296  double v = m_pxdclusters[iclus]->getV();
297  TVector3 current(u, v, 0);
298 
299  //2D dist sqared
300  double dist = (intersection - current).Mag();
301  if (dist < mindist) {
302  closest = iclus;
303  mindist = dist;
304  }
305  }
306 
307  return closest;
308 
309 }
310 
311 bool PXDDQMEfficiencyModule::isCloseToBorder(int u, int v, int checkDistance)
312 {
313 
314  if (u - checkDistance < 0 || u + checkDistance >= 250 ||
315  v - checkDistance < 0 || v + checkDistance >= 768) {
316  return true;
317  }
318  return false;
319 }
320 
321 bool PXDDQMEfficiencyModule::isDeadPixelClose(int u, int v, int checkDistance, const VxdID& moduleID)
322 {
323 
324  //Iterate over square around the intersection to see if any close pixel is dead
325  for (int u_iter = u - checkDistance; u_iter <= u + checkDistance ; ++u_iter) {
326  for (int v_iter = v - checkDistance; v_iter <= v + checkDistance ; ++v_iter) {
327  if (PXD::PXDPixelMasker::getInstance().pixelDead(moduleID, u_iter, v_iter)
328  || !PXD::PXDPixelMasker::getInstance().pixelOK(moduleID, u_iter, v_iter)) {
329  return true;
330  }
331  }
332  }
333  return false;
334 }
Belle2::RelationVector::size
size_t size() const
Get number of relations.
Definition: RelationVector.h:98
Belle2::PXDDQMEfficiencyModule::event
void event() override final
main function which fills trees and histograms
Definition: PXDDQMEfficiencyModule.cc:95
Belle2::PXDDQMEfficiencyModule::m_vxdGeometry
VXD::GeoCache & m_vxdGeometry
the geometry
Definition: PXDDQMEfficiencyModule.h:109
Belle2::PXDDQMEfficiencyModule::m_h_matched_cluster
std::map< VxdID, TH2F * > m_h_matched_cluster
histograms of matched clusters
Definition: PXDDQMEfficiencyModule.h:141
Belle2::PXDDQMEfficiencyModule::isCloseToBorder
bool isCloseToBorder(int u, int v, int checkDistance)
is it close to the border
Definition: PXDDQMEfficiencyModule.cc:311
Belle2::PXDDQMEfficiencyModule::m_z0maxCut
double m_z0maxCut
cut z0 maximum in cm (large positiv value eg 9999 disables)
Definition: PXDDQMEfficiencyModule.h:135
Belle2::VxdID
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:43
Belle2::PXDDQMEfficiencyModule::m_histogramDirectoryName
std::string m_histogramDirectoryName
Where to save the histograms too.
Definition: PXDDQMEfficiencyModule.h:112
genfit::MeasuredStateOnPlane
#StateOnPlane with additional covariance matrix.
Definition: MeasuredStateOnPlane.h:39
Belle2::PXDDQMEfficiencyModule::m_h_pt2
std::map< VxdID, TH1F * > m_h_pt2
histograms of pt2
Definition: PXDDQMEfficiencyModule.h:147
Belle2::PXDDQMEfficiencyModule::m_maskDeadPixels
bool m_maskDeadPixels
mask dead pixels
Definition: PXDDQMEfficiencyModule.h:102
Belle2::PXDDQMEfficiencyModule::m_momCut
double m_momCut
Cut on fitted track momentum.
Definition: PXDDQMEfficiencyModule.h:131
Belle2::PXDDQMEfficiencyModule::m_PXDInterceptListName
std::string m_PXDInterceptListName
intercept list name
Definition: PXDDQMEfficiencyModule.h:117
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::PXDDQMEfficiencyModule::m_requireROIs
bool m_requireROIs
Require tracks going through ROIs.
Definition: PXDDQMEfficiencyModule.h:97
Belle2::PXDDQMEfficiencyModule::m_z0minCut
double m_z0minCut
cut z0 minimum in cm (large negativ value eg -9999 disables)
Definition: PXDDQMEfficiencyModule.h:134
Belle2::VxdID::getLadderNumber
baseType getLadderNumber() const
Get the ladder id.
Definition: VxdID.h:108
Belle2::VXD::SensorInfoBase
Base class to provide Sensor Information for PXD and SVD.
Definition: SensorInfoBase.h:40
Belle2::PXDDQMEfficiencyModule::beginRun
void beginRun() override final
begin run function which resets histograms
Definition: PXDDQMEfficiencyModule.cc:81
Belle2::PXDDQMEfficiencyModule::m_tracks
StoreArray< RecoTrack > m_tracks
store array of tracks
Definition: PXDDQMEfficiencyModule.h:123
Belle2::PXDDQMEfficiencyModule::m_h_sv
std::map< VxdID, TH1F * > m_h_sv
histograms of sv
Definition: PXDDQMEfficiencyModule.h:145
Belle2::VXD::GeoCache::getListOfSensors
const std::vector< VxdID > getListOfSensors() const
Get list of all sensors.
Definition: GeoCache.cc:60
Belle2::PXDDQMEfficiencyModule::m_distcut
double m_distcut
distance cut in cm!
Definition: PXDDQMEfficiencyModule.h:127
Belle2::PXDDQMEfficiencyModule::m_h_track_hits
std::map< VxdID, TH2F * > m_h_track_hits
histograms of track hits
Definition: PXDDQMEfficiencyModule.h:140
Belle2::PXDDQMEfficiencyModule::m_verboseHistos
bool m_verboseHistos
add some verbose histograms for cuts
Definition: PXDDQMEfficiencyModule.h:106
Belle2::PXDDQMEfficiencyModule::m_tracksName
std::string m_tracksName
name of the store array of tracks
Definition: PXDDQMEfficiencyModule.h:115
Belle2::PXDDQMEfficiencyModule::m_pcut
double m_pcut
pValue-Cut for tracks
Definition: PXDDQMEfficiencyModule.h:130
Belle2::PXDDQMEfficiencyModule::m_ROIs
StoreArray< ROIid > m_ROIs
store array of ROIs
Definition: PXDDQMEfficiencyModule.h:124
Belle2::Const::pion
static const ChargedStable pion
charged pion particle
Definition: Const.h:535
genfit::FitStatus::getPVal
virtual double getPVal() const
Get the p value of the fit.
Definition: FitStatus.h:128
Belle2::PXDDQMEfficiencyModule::m_h_pt
std::map< VxdID, TH1F * > m_h_pt
histograms of transverse momenta
Definition: PXDDQMEfficiencyModule.h:143
Belle2::RelationVector
Class for type safe access to objects that are referred to in relations.
Definition: DataStore.h:38
Belle2::PXDDQMEfficiencyModule::m_ROIsName
std::string m_ROIsName
name of the store array of ROIs
Definition: PXDDQMEfficiencyModule.h:116
Belle2::PXDDQMEfficiencyModule::m_h_p
std::map< VxdID, TH1F * > m_h_p
histograms of momenta
Definition: PXDDQMEfficiencyModule.h:142
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::PXDDQMEfficiencyModule::m_pxdclusters
StoreArray< PXDCluster > m_pxdclusters
store array of pxd clusters
Definition: PXDDQMEfficiencyModule.h:122
Belle2::PXDDQMEfficiencyModule::m_cutBorders
bool m_cutBorders
cut borders
Definition: PXDDQMEfficiencyModule.h:104
Belle2::PXDDQMEfficiencyModule::m_h_p2
std::map< VxdID, TH1F * > m_h_p2
histograms of p2
Definition: PXDDQMEfficiencyModule.h:146
Belle2::PXDDQMEfficiencyModule::m_vFactor
double m_vFactor
factor for track-error on distcut comparison
Definition: PXDDQMEfficiencyModule.h:129
Belle2::PXDDQMEfficiencyModule::isDeadPixelClose
bool isDeadPixelClose(int u, int v, int checkDistance, const VxdID &moduleID)
is a dead pixel close
Definition: PXDDQMEfficiencyModule.cc:321
Belle2::PXDDQMEfficiencyModule::m_uFactor
double m_uFactor
factor for track-error on distcut comparison
Definition: PXDDQMEfficiencyModule.h:128
Belle2::PXDDQMEfficiencyModule::findClosestCluster
int findClosestCluster(const VxdID &vxdid, TVector3 intersection)
find the closest cluster
Definition: PXDDQMEfficiencyModule.cc:276
Belle2::PXDDQMEfficiencyModule::m_pxdClustersName
std::string m_pxdClustersName
name of the store array of pxd clusters
Definition: PXDDQMEfficiencyModule.h:114
Belle2::PXDDQMEfficiencyModule::m_u_bins
int m_u_bins
the u bins
Definition: PXDDQMEfficiencyModule.h:119
Belle2::PXDDQMEfficiencyModule::m_intercepts
StoreArray< PXDIntercept > m_intercepts
store array of PXD Intercepts
Definition: PXDDQMEfficiencyModule.h:125
Belle2::PXDDQMEfficiencyModule::initialize
void initialize() override final
initializes the need store arrays, trees and histograms
Definition: PXDDQMEfficiencyModule.cc:68
Belle2::StoreArray::isValid
bool isValid() const
Check wether the array was registered.
Definition: StoreArray.h:298
Belle2::VxdID::getSensorNumber
baseType getSensorNumber() const
Get the sensor id.
Definition: VxdID.h:110
Belle2::VXD::SensorInfoBase::PXD
@ PXD
PXD Sensor.
Definition: SensorInfoBase.h:44
Belle2::PXDDQMEfficiencyModule::m_d0Cut
double m_d0Cut
cut abs(d0) in cm (and negativ value eg -9999 disables)
Definition: PXDDQMEfficiencyModule.h:136
Belle2::PXDDQMEfficiencyModule::m_h_sv2
std::map< VxdID, TH1F * > m_h_sv2
histrograms of sv2
Definition: PXDDQMEfficiencyModule.h:149
Belle2::PXDDQMEfficiencyModule::defineHisto
void defineHisto() override final
actually defines the trees and histograms
Definition: PXDDQMEfficiencyModule.cc:232
Belle2::PXDDQMEfficiencyModule::m_h_su
std::map< VxdID, TH1F * > m_h_su
histograms of su
Definition: PXDDQMEfficiencyModule.h:144
Belle2::PXDDQMEfficiencyModule::m_v_bins
int m_v_bins
the v bins
Definition: PXDDQMEfficiencyModule.h:120
Belle2::PXDDQMEfficiencyModule::m_pTCut
double m_pTCut
Cut on fitted track pT.
Definition: PXDDQMEfficiencyModule.h:132
Belle2::VXD::GeoCache::getSensorInfo
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:68
Belle2::Track
Class that bundles various TrackFitResults.
Definition: Track.h:35
Belle2::VxdID::getLayerNumber
baseType getLayerNumber() const
Get the layer id.
Definition: VxdID.h:106
Belle2::PXDDQMEfficiencyModule
Creates the basic histograms for PXD Efficiency DQM Simplified and adopted version of the testbeam px...
Definition: PXDDQMEfficiencyModule.h:54
genfit::FitStatus
Class where important numbers and properties of a fit can be stored.
Definition: FitStatus.h:80
Belle2::PXDIntercept
PXDIntercept stores the U,V coordinates and uncertainties of the intersection of a track with an PXD ...
Definition: PXDIntercept.h:32
Belle2::PXDDQMEfficiencyModule::m_maskedDistance
int m_maskedDistance
Distance inside which no dead pixel or module border is allowed.
Definition: PXDDQMEfficiencyModule.h:137
Belle2::PXDDQMEfficiencyModule::m_h_su2
std::map< VxdID, TH1F * > m_h_su2
histrograms of su2
Definition: PXDDQMEfficiencyModule.h:148
Belle2::StoreArray::getEntries
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:226
Belle2::HistoModule
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
Belle2::PXD::PXDPixelMasker::getInstance
static PXDPixelMasker & getInstance()
Main (and only) way to access the PXDPixelMasker.
Definition: PXDPixelMasker.cc:37
Belle2::PXDDQMEfficiencyModule::m_minSVDHits
unsigned int m_minSVDHits
Required hits in SVD strips for tracks.
Definition: PXDDQMEfficiencyModule.h:133