Belle II Software  release-06-02-00
PXDDQMEfficiencySelftrackModule.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/PXDDQMEfficiencySelftrackModule.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(PXDDQMEfficiencySelftrack)
24 
25 //-----------------------------------------------------------------
26 // Implementation
27 //-----------------------------------------------------------------
28 
29 PXDDQMEfficiencySelftrackModule::PXDDQMEfficiencySelftrackModule() : 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?
35  setPropertyFlags(c_ParallelProcessingCertified);
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("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms will be placed",
42  std::string("PXDEFF"));
43  addParam("binsU", m_u_bins, "histogram bins in u direction", int(4));
44  addParam("binsV", m_v_bins, "histogram bins in v direction", int(6));
45  addParam("distCut", m_distcut, "max distance in [cm] for cluster to be counted to a track", double(0.0500));
46  addParam("pCut", m_pcut, "Set a cut on the track p-value", double(1e-20));
47  addParam("requireROIs", m_requireROIs, "require tracks to lie inside a ROI", bool(false));
48  addParam("useAlignment", m_useAlignment, "if true the alignment will be used", bool(true));
49  addParam("maskDeadPixels", m_maskDeadPixels, "Do not consider tracks going through known dead or hot pixels for the efficiency",
50  bool(false));
51  addParam("minSVDHits", m_minSVDHits, "Number of SVD hits required in a track to be considered", 0u);
52  addParam("momCut", m_momCut, "Set a cut on the track momentum in GeV/c, 0 disables", double(0));
53  addParam("pTCut", m_pTCut, "Set a cut on the track pT in GeV/c, 0 disables", double(1));
54  addParam("cutBorders", m_cutBorders, "Do not use tracks near the borders of the sensor", bool(true));
55  addParam("maskedDistance", m_maskedDistance, "Distance inside which no masked pixel or sensor border is allowed", int(10));
56  addParam("trackUFactorDistCut", m_uFactor, "Set a cut on u error of track (factor*err<dist), 0 disables", double(2.0));
57  addParam("trackVFactorDistCut", m_vFactor, "Set a cut on v error of track (factor*err<dist), 0 disables", double(2.0));
58  addParam("z0minCut", m_z0minCut, "Set a cut z0 minimum in cm (large negativ value eg -9999 disables)", double(-1));
59  addParam("z0maxCut", m_z0maxCut, "Set a cut z0 maximum in cm (large positiv value eg 9999 disables)", double(1));
60  addParam("d0Cut", m_d0Cut, "Set a cut abs(d0) in cm (and negativ value eg -9999 disables)", double(0.5));
61  addParam("verboseHistos", m_verboseHistos, "Add more verbose histograms for cuts (not for ereoc)", bool(false));
62 }
63 
64 
66 {
67  //calls the define histogram function
68  REG_HISTOGRAM;
69 
70  //register the required arrays
71  //Register as optional so validation for cases where they are not available still succeeds, but module will not do any meaningful work without them
73  m_tracks.isOptional(m_tracksName);
74  m_ROIs.isOptional(m_ROIsName);
75 }
76 
78 {
79  for (auto& h : m_h_track_hits) if (h.second) h.second->Reset();
80  for (auto& h : m_h_matched_cluster) if (h.second) h.second->Reset();
81  for (auto& h : m_h_p) if (h.second) h.second->Reset();
82  for (auto& h : m_h_pt) if (h.second) h.second->Reset();
83  for (auto& h : m_h_su) if (h.second) h.second->Reset();
84  for (auto& h : m_h_sv) if (h.second) h.second->Reset();
85  for (auto& h : m_h_p2) if (h.second) h.second->Reset();
86  for (auto& h : m_h_pt2) if (h.second) h.second->Reset();
87  for (auto& h : m_h_su2) if (h.second) h.second->Reset();
88  for (auto& h : m_h_sv2) if (h.second) h.second->Reset();
89 }
90 
92 {
93  if (!m_pxdclusters.isValid()) {
94  B2INFO("PXDClusters array is missing, no efficiencies");
95  return;
96  }
97  if (!m_tracks.isValid()) {
98  B2INFO("RecoTrack array is missing, no efficiencies");
99  return;
100  }
101  if (!m_ROIs.isValid() && m_requireROIs) {
102  B2INFO("ROI array is missing but required hits in ROIs, aborting");
103  return;
104  }
105 
106 
107  for (auto& a_track : m_tracks) {
108 
109  //If fit failed assume position pointed to is useless anyway
110  if (!a_track.wasFitSuccessful()) continue;
111 
112  if (a_track.getNumberOfSVDHits() < m_minSVDHits) continue;
113 
114  const genfit::FitStatus* fitstatus = a_track.getTrackFitStatus();
115  if (fitstatus->getPVal() < m_pcut) continue;
116 
117  genfit::MeasuredStateOnPlane trackstate;
118  trackstate = a_track.getMeasuredStateOnPlaneFromFirstHit();
119  if (trackstate.getMom().Mag() < m_momCut) continue;
120  if (trackstate.getMom().Pt() < m_pTCut) continue;
121 
122  auto ptr = a_track.getRelated<Track>("Tracks");
123 
124  if (!ptr) {
125  B2ERROR("expect a track for fitted recotracks");
126  continue;
127  }
128  auto ptr2 = ptr->getTrackFitResultWithClosestMass(Const::pion);
129  if (!ptr2) {
130  B2ERROR("expect a track fit result for mass");
131  continue;
132  }
133 
134  // Vertex cut
135  if (ptr2->getZ0() < m_z0minCut || ptr2->getZ0() > m_z0maxCut || fabs(ptr2->getD0()) > m_d0Cut) continue;
136 
137  //loop over all PXD sensors to get the intersections
138  std::vector<VxdID> sensors = m_vxdGeometry.getListOfSensors();
139  for (VxdID& aVxdID : sensors) {
141  if (info.getType() != VXD::SensorInfoBase::PXD) continue;
142  //Search for intersections of the track with all PXD layers
143  //Traditional (aka the person before did it like this) method
144  //If there is a way to find out sensors crossed by a track directly, that would most likely be faster
145 
146  bool isgood = false;
147  //true = track intersects current sensor
148  double sigu(-9999);
149  double sigv(-9999);
150  TVector3 intersec_buff = getTrackInterSec(info, a_track, isgood, sigu, sigv);
151 
152  if (!isgood) {
153  continue;//track does not go through this sensor-> nothing to measure anyway
154  } else {
155  if (m_verboseHistos) {
156  if (m_h_p[aVxdID]) m_h_p[aVxdID]->Fill(trackstate.getMom().Mag());
157  if (m_h_pt[aVxdID]) m_h_pt[aVxdID]->Fill(trackstate.getMom().Pt());
158  if (m_h_su[aVxdID]) m_h_su[aVxdID]->Fill(sigu);
159  if (m_h_sv[aVxdID]) m_h_sv[aVxdID]->Fill(sigv);
160  }
161  if (m_uFactor * sigu > m_distcut) continue; // Error ufak*SigmaU > cut
162  if (m_vFactor * sigv > m_distcut) continue; // Error vfak*SigmaV > cut
163 
164  double u_fit = intersec_buff.X();
165  double v_fit = intersec_buff.Y();
166 
167  int ucell_fit = info.getUCellID(intersec_buff.X());
168  int vcell_fit = info.getVCellID(intersec_buff.Y());
169 
170  if (m_cutBorders && isCloseToBorder(ucell_fit, vcell_fit, m_maskedDistance)) {
171  continue;
172  }
173 
174  if (m_maskDeadPixels && isDeadPixelClose(ucell_fit, vcell_fit, m_maskedDistance, aVxdID)) {
175  continue;
176  }
177 
178  if (m_requireROIs) {
179  //Check if the intersection is inside a ROI
180  //If not, even if measured the cluster was thrown away->Not PXD's fault
181  bool fitInsideROI = false;
182  for (auto& roit : m_ROIs) {
183  if (aVxdID != roit.getSensorID()) {
184  continue; //ROI on other sensor
185  }
186 
187  if (ucell_fit < roit.getMaxUid()
188  && ucell_fit > roit.getMinUid()
189  && vcell_fit < roit.getMaxVid()
190  && vcell_fit > roit.getMinVid()) {
191  fitInsideROI = true;
192  }
193  }
194  if (!fitInsideROI) {
195  continue;//Hit wouldn't have been recorded
196  }
197  }
198 
199  //This track should be on the sensor
200  m_h_track_hits[aVxdID]->Fill(ucell_fit, vcell_fit);
201 
202  //Now check if the sensor measured a hit here
203 
204  int bestcluster = findClosestCluster(aVxdID, intersec_buff);
205  if (bestcluster >= 0) {
206  double u_clus = m_pxdclusters[bestcluster]->getU();
207  double v_clus = m_pxdclusters[bestcluster]->getV();
208 
209  //is the closest cluster close enough to the track to count as measured?
210  TVector3 dist_clus(u_fit - u_clus, v_fit - v_clus, 0);
211  if (dist_clus.Mag() <= m_distcut) {
212  m_h_matched_cluster[aVxdID]->Fill(ucell_fit, vcell_fit);
213  if (m_verboseHistos) {
214  if (m_h_p2[aVxdID]) m_h_p2[aVxdID]->Fill(trackstate.getMom().Mag());
215  if (m_h_pt2[aVxdID]) m_h_pt2[aVxdID]->Fill(trackstate.getMom().Pt());
216  if (m_h_su2[aVxdID]) m_h_su2[aVxdID]->Fill(sigu);
217  if (m_h_sv2[aVxdID]) m_h_sv2[aVxdID]->Fill(sigv);
218  }
219  }
220  }
221  }
222  }
223  }
224 }
225 
226 
227 
229  bool& isgood,
230  double& du, double& dv)
231 {
232  //will be set true if the intersect was found
233  isgood = false;
234 
235  TVector3 intersec(99999999, 9999999, 0); //point outside the sensor
236 
238 
239  //adopted (aka stolen) from tracking/modules/pxdClusterRescue/PXDClusterRescueROIModule
240  try {
241  // get sensor plane
242  TVector3 zeroVec(0, 0, 0);
243  TVector3 uVec(1, 0, 0);
244  TVector3 vVec(0, 1, 0);
245 
246  genfit::DetPlane* sensorPlane = new genfit::DetPlane();
247  sensorPlane->setO(pxdSensorInfo.pointToGlobal(zeroVec, m_useAlignment));
248  sensorPlane->setUV(pxdSensorInfo.vectorToGlobal(uVec, m_useAlignment), pxdSensorInfo.vectorToGlobal(vVec, m_useAlignment));
249 
250  //boost pointer (will be deleted automatically ?!?!?)
251  genfit::SharedPlanePtr sensorPlaneSptr(sensorPlane);
252 
253  // do extrapolation
254  gfTrackState.extrapolateToPlane(sensorPlaneSptr);
255  } catch (genfit::Exception& gfException) {
256  B2WARNING("Fitting failed: " << gfException.getExcString());
257  isgood = false;
258  return intersec;
259  } catch (...) {
260  B2WARNING("Fitting failed: for some reason");
261  isgood = false;
262  return intersec;
263  }
264 
265  //local position
266  intersec = pxdSensorInfo.pointToLocal(gfTrackState.getPos(), m_useAlignment);
267 
268  //try to get the momentum
269  B2DEBUG(1, "Fitted momentum on the plane p = " << gfTrackState.getMom().Mag());
270 
271  // no tolerance currently! Maybe one should be added!
272  double tolerance = 0.0;
273  bool inside = pxdSensorInfo.inside(intersec.X(), intersec.Y(), tolerance, tolerance);
274 
275  // get intersection point in local coordinates with covariance matrix
276  TMatrixDSym covMatrix = gfTrackState.getCov(); // 5D with elements q/p,u',v',u,v in plane system
277 
278  // get ROI by covariance matrix and local intersection point
279  du = std::sqrt(covMatrix(3, 3));
280  dv = std::sqrt(covMatrix(4, 4));
281 
282  if (inside) isgood = true;
283 
284  return intersec;
285 }
286 
287 
289 {
290  // Create a separate histogram directory and cd into it.
291  TDirectory* oldDir = gDirectory;
292  if (m_histogramDirectoryName != "") {
293  oldDir->mkdir(m_histogramDirectoryName.c_str());
294  oldDir->cd(m_histogramDirectoryName.c_str());
295  }
296 
297  std::vector<VxdID> sensors = m_vxdGeometry.getListOfSensors();
298  for (VxdID& avxdid : sensors) {
300  if (info.getType() != VXD::SensorInfoBase::PXD) continue;
301  //Only interested in PXD sensors
302 
303  TString buff = (std::string)avxdid;
304  buff.ReplaceAll(".", "_");
305 
306  int nu = info.getUCells();
307  int nv = info.getVCells();
308 
309  //nu + 1,nv + 1 Bin numbers when using one bin per pixel
310  m_h_track_hits[avxdid] = new TH2F("track_hits_" + buff, "tracks through sensor " + buff,
311  m_u_bins, -0.5, nu - 0.5, m_v_bins, -0.5, nv - 0.5);
312  m_h_matched_cluster[avxdid] = new TH2F("matched_cluster_" + buff, "clusters matched to track intersections " + buff,
313  m_u_bins, -0.5, nu - 0.5, m_v_bins, -0.5, nv - 0.5);
314 
315  if (m_verboseHistos) {
316  m_h_p[avxdid] = new TH1F("p_" + buff, "p " + buff, 100, 0, 10);
317  m_h_pt[avxdid] = new TH1F("pt_" + buff, "pt " + buff, 100, 0, 10);
318  m_h_su[avxdid] = new TH1F("su_" + buff, "su " + buff, 1000, 0, 1);
319  m_h_sv[avxdid] = new TH1F("sv_" + buff, "sv " + buff, 1000, 0, 1);
320  m_h_p2[avxdid] = new TH1F("p2_" + buff, "p2 " + buff, 100, 0, 10);
321  m_h_pt2[avxdid] = new TH1F("pt2_" + buff, "pt2 " + buff, 100, 0, 10);
322  m_h_su2[avxdid] = new TH1F("su2_" + buff, "su2 " + buff, 1000, 0, 1);
323  m_h_sv2[avxdid] = new TH1F("sv2_" + buff, "sv2 " + buff, 1000, 0, 1);
324  }
325  }
326  // cd back to root directory
327  oldDir->cd();
328 }
329 
330 
331 int
332 PXDDQMEfficiencySelftrackModule::findClosestCluster(const VxdID& avxdid, TVector3 intersection)
333 {
334  int closest = -1;
335  double mindist = 999999999999; //definitely outside of the sensor
336 
338 
339  //loop the clusters
340  for (int iclus = 0; iclus < m_pxdclusters.getEntries(); iclus++) {
341  //Do not consider as different if only segment differs!
342  //As of this writing segment is never filled for clusters, but just to be sure
343  VxdID clusterID = m_pxdclusters[iclus]->getSensorID();
344  if (avxdid.getLayerNumber() != clusterID.getLayerNumber() ||
345  avxdid.getLadderNumber() != clusterID.getLadderNumber() ||
346  avxdid.getSensorNumber() != clusterID.getSensorNumber()) {
347  continue;
348  }
349  //only cluster on the correct sensor and direction should survive
350 
351  double u = m_pxdclusters[iclus]->getU();
352  double v = m_pxdclusters[iclus]->getV();
353  TVector3 current(u, v, 0);
354 
355  //2D dist sqared
356  double dist = (intersection - current).Mag();
357  if (dist < mindist) {
358  closest = iclus;
359  mindist = dist;
360  }
361  }
362 
363  return closest;
364 
365 }
366 
367 bool PXDDQMEfficiencySelftrackModule::isCloseToBorder(int u, int v, int checkDistance)
368 {
369 
370  if (u - checkDistance < 0 || u + checkDistance >= 250 ||
371  v - checkDistance < 0 || v + checkDistance >= 768) {
372  return true;
373  }
374  return false;
375 }
376 
377 bool PXDDQMEfficiencySelftrackModule::isDeadPixelClose(int u, int v, int checkDistance, const VxdID& moduleID)
378 {
379 
380  //Iterate over square around the intersection to see if any close pixel is dead
381  for (int u_iter = u - checkDistance; u_iter <= u + checkDistance ; ++u_iter) {
382  for (int v_iter = v - checkDistance; v_iter <= v + checkDistance ; ++v_iter) {
383  if (PXD::PXDPixelMasker::getInstance().pixelDead(moduleID, u_iter, v_iter)
384  || !PXD::PXDPixelMasker::getInstance().pixelOK(moduleID, u_iter, v_iter)) {
385  return true;
386  }
387  }
388  }
389  return false;
390 }
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
Creates the basic histograms for PXD Efficiency DQM Simplified and adopted version of the testbeam px...
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 (and negativ value eg -9999 disables)
TVector3 getTrackInterSec(const VXD::SensorInfoBase &pxdSensorInfo, const RecoTrack &aTrack, bool &isgood, double &du, double &dv)
helper functions to do some of the calculations
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
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::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)
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.
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.
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
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
static PXDPixelMasker & getInstance()
Main (and only) way to access the PXDPixelMasker.
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:76
const genfit::MeasuredStateOnPlane & getMeasuredStateOnPlaneFromFirstHit(const genfit::AbsTrackRep *representation=nullptr) const
Return genfit's MeasuredStateOnPlane for the first hit in a fit useful for extrapolation of measureme...
Definition: RecoTrack.cc:586
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:58
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:66
Base class to provide Sensor Information for PXD and SVD.
bool inside(double u, double v, double uTolerance=DBL_EPSILON, double vTolerance=DBL_EPSILON) const
Check wether a given point is inside the active area.
TVector3 pointToLocal(const TVector3 &global, bool reco=false) const
Convert a point from global to local coordinates.
TVector3 pointToGlobal(const TVector3 &local, bool reco=false) const
Convert a point from local to global coordinates.
TVector3 vectorToGlobal(const TVector3 &local, bool reco=false) const
Convert a vector from local to global coordinates.
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
Detector plane.
Definition: DetPlane.h:59
Exception class for error handling in GENFIT (provides storage for diagnostic information)
Definition: Exception.h:48
Class where important numbers and properties of a fit can be stored.
Definition: FitStatus.h:80
virtual double getPVal() const
Get the p value of the fit.
Definition: FitStatus.h:128
#StateOnPlane with additional covariance matrix.
#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.
std::shared_ptr< genfit::DetPlane > SharedPlanePtr
Shared Pointer to a DetPlane.