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