Belle II Software  release-06-02-00
PXDDQMEfficiencyNtupleSelftrackModule.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/PXDDQMEfficiencyNtupleSelftrackModule.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 "TMatrixDSym.h"
17 using namespace Belle2;
18 
19 //-----------------------------------------------------------------
20 // Register the Module
21 //-----------------------------------------------------------------
22 REG_MODULE(PXDDQMEfficiencyNtupleSelftrack)
23 
24 //-----------------------------------------------------------------
25 // Implementation
26 //-----------------------------------------------------------------
27 
29  m_vxdGeometry(VXD::GeoCache::getInstance())
30 {
31  // Set module properties
32  setDescription("Create basic histograms for PXD efficiency");
33 
34  // setPropertyFlags(c_ParallelProcessingCertified);// for ntuple not certified!!!
35 
36  // Parameter definitions
37  addParam("pxdClustersName", m_pxdClustersName, "name of StoreArray with PXD cluster", std::string(""));
38  addParam("recoTracksName", m_recoTracksName, "name of StoreArray with RecoTracks", std::string(""));
39  addParam("tracksName", m_tracksName, "name of StoreArray with Tracks", std::string(""));
40  addParam("ROIsName", m_ROIsName, "name of the list of HLT ROIs, if available in output", std::string(""));
41  addParam("useAlignment", m_useAlignment, "if true the alignment will be used", true);
42  addParam("pCut", m_pcut, "Set a cut on the track fit p-value (0=no cut)", double(0));
43  addParam("minSVDHits", m_minSVDHits, "Number of SVD hits required in a track to be considered", 0u);
44  addParam("momCut", m_momCut, "Set a cut on the track momentum in GeV/c, 0 disables", double(0));
45  addParam("pTCut", m_pTCut, "Set a cut on the track pT in GeV/c, 0 disables", double(0));
46  addParam("maskedDistance", m_maskedDistance, "Distance inside which no masked pixel or sensor border is allowed", int(10));
47 }
48 
49 
51 {
52  auto dir = gDirectory;
53  if (m_tuple) {
54  if (m_file) { // no file -> no write to file
55  m_file->cd();
56  }
57  m_tuple->Write();
58  delete m_tuple;
59  m_tuple = nullptr;
60  }
61  if (m_file) {
62  m_file->Write();
63  m_file->Close();
64  delete m_file;
65  m_file = nullptr;
66  }
67  dir->cd();
68 }
69 
70 
72 {
73  m_file = new TFile("test.root", "recreate");
74  if (m_file) m_file->cd();
75  m_tuple = new TNtuple("effcontrol", "effcontrol",
76  "vxdid:u:v:p:pt:distu:distv:sigu:sigv:dist:inroi:clborder:cldead:matched:z0:d0:svdhits:charge:phi:costheta");
77 
78  //register the required arrays
79  //Register as optional so validation for cases where they are not available still succeeds, but module will not do any meaningful work without them
81  m_recoTracks.isOptional(m_recoTracksName);
82  m_tracks.isOptional(m_tracksName);
83  m_ROIs.isOptional(m_ROIsName);
84 }
85 
86 
88 {
89  if (!m_pxdclusters.isValid()) {
90  B2INFO("PXDClusters array is missing, no efficiencies");
91  return;
92  }
93  if (!m_recoTracks.isValid()) {
94  B2INFO("RecoTrack array is missing, no efficiencies");
95  return;
96  }
97  if (!m_tracks.isValid()) {
98  B2INFO("Track array is missing, no efficiencies");
99  return;
100  }
101 
102  for (auto& track : m_tracks) {
103  RelationVector<RecoTrack> recoTrack = track.getRelationsTo<RecoTrack>(m_recoTracksName);
104  if (!recoTrack.size()) continue;
105 
106  auto a_track = recoTrack[0];
107  //If fit failed assume position pointed to is useless anyway
108  if (!a_track->wasFitSuccessful()) continue;
109 
110  if (a_track->getNumberOfSVDHits() < m_minSVDHits) continue;
111 
112  const genfit::FitStatus* fitstatus = a_track->getTrackFitStatus();
113  if (fitstatus->getPVal() < m_pcut) continue;
114 
115  genfit::MeasuredStateOnPlane trackstate;
116  trackstate = a_track->getMeasuredStateOnPlaneFromFirstHit();
117  if (trackstate.getMom().Mag() < m_momCut) continue;
118  if (trackstate.getMom().Pt() < m_pTCut) continue;
119 
120  const TrackFitResult* ptr2 = track.getTrackFitResultWithClosestMass(Const::pion);
121  if (!ptr2) {
122  B2ERROR("expect a track fit result for mass");
123  continue;
124  }
125 
126  //loop over all PXD sensors to get the intersections
127  std::vector<VxdID> sensors = m_vxdGeometry.getListOfSensors();
128  for (VxdID& aVxdID : sensors) {
130  if (info.getType() != VXD::SensorInfoBase::PXD) continue;
131  //Search for intersections of the track with all PXD layers
132  //Traditional (aka the person before did it like this) method
133  //If there is a way to find out sensors crossed by a track directly, that would most likely be faster
134 
135  bool isgood = false;
136  //true = track intersects current sensor
137  double sigu(-9999);
138  double sigv(-9999);
139  TVector3 intersec_buff = getTrackInterSec(info, *a_track, isgood, sigu, sigv);
140 
141  if (!isgood) {
142  continue;//track does not go through this sensor-> nothing to measure anyway
143  } else {
144 
145  double u_fit = intersec_buff.X();
146  double v_fit = intersec_buff.Y();
147 
148  int ucell_fit = info.getUCellID(intersec_buff.X());
149  int vcell_fit = info.getVCellID(intersec_buff.Y());
150 
151  bool closeToBoarder = false;
152  if (isCloseToBorder(ucell_fit, vcell_fit, m_maskedDistance)) {
153  closeToBoarder = true;
154  }
155 
156  bool closeToDead = false;
157  if (isDeadPixelClose(ucell_fit, vcell_fit, m_maskedDistance, aVxdID)) {
158  closeToDead = true;
159  }
160 
161  bool fitInsideROI = false;
162  if (m_ROIs.isValid()) {
163  //Check if the intersection is inside a ROI
164  //If not, even if measured the cluster was thrown away->Not PXD's fault
165  for (auto& roit : m_ROIs) {
166  if (aVxdID != roit.getSensorID()) {
167  continue; //ROI on other sensor
168  }
169 
170  if (ucell_fit < roit.getMaxUid()
171  && ucell_fit > roit.getMinUid()
172  && vcell_fit < roit.getMaxVid()
173  && vcell_fit > roit.getMinVid()) {
174  fitInsideROI = true;
175  }
176  }
177  }
178 
179  //Now check if the sensor measured a hit here
180 
181  int bestcluster = findClosestCluster(aVxdID, intersec_buff);
182  double du_clus = 0;
183  double dv_clus = 0;
184  double d_clus = 0;
185  float charge = 0;
186  bool matched = false;
187  if (bestcluster >= 0) {
188  double u_clus = m_pxdclusters[bestcluster]->getU();
189  double v_clus = m_pxdclusters[bestcluster]->getV();
190 
191  //is the closest cluster close enough to the track to count as measured?
192  TVector3 dist_clus(u_fit - u_clus, v_fit - v_clus, 0);
193  du_clus = u_fit - u_clus;
194  dv_clus = v_fit - v_clus;
195  d_clus = dist_clus.Mag();
196  charge = m_pxdclusters[bestcluster]->getCharge();
197  matched = true;
198  }
199  float fill[22] = {float((int)aVxdID), float(u_fit), float(v_fit), float(trackstate.getMom().Mag()), float(trackstate.getMom().Pt()),
200  float(du_clus), float(dv_clus), float(sigu), float(sigv), float(d_clus),
201  float(fitInsideROI), float(closeToBoarder), float(closeToDead), float(matched),
202  float(ptr2->getZ0()), float(ptr2->getD0()), float(a_track->getNumberOfSVDHits()),
203  charge, float(trackstate.getMom().Phi()), float(trackstate.getMom().CosTheta())
204  };
205  m_tuple->Fill(fill);
206  }
207  }
208  }
209 }
210 
211 
212 
214  bool& isgood,
215  double& du, double& dv)
216 {
217  //will be set true if the intersect was found
218  isgood = false;
219 
220  TVector3 intersec(99999999, 9999999, 0); //point outside the sensor
221 
223 
224  //adopted (aka stolen) from tracking/modules/pxdClusterRescue/PXDClusterRescueROIModule
225  try {
226  // get sensor plane
227  TVector3 zeroVec(0, 0, 0);
228  TVector3 uVec(1, 0, 0);
229  TVector3 vVec(0, 1, 0);
230 
231  genfit::DetPlane* sensorPlane = new genfit::DetPlane();
232  sensorPlane->setO(pxdSensorInfo.pointToGlobal(zeroVec, m_useAlignment));
233  sensorPlane->setUV(pxdSensorInfo.vectorToGlobal(uVec, m_useAlignment), pxdSensorInfo.vectorToGlobal(vVec, m_useAlignment));
234 
235  //boost pointer (will be deleted automatically ?!?!?)
236  genfit::SharedPlanePtr sensorPlaneSptr(sensorPlane);
237 
238  // do extrapolation
239  gfTrackState.extrapolateToPlane(sensorPlaneSptr);
240  } catch (genfit::Exception& gfException) {
241  B2WARNING("Fitting failed: " << gfException.getExcString());
242  isgood = false;
243  return intersec;
244  } catch (...) {
245  B2WARNING("Fitting failed: for some reason");
246  isgood = false;
247  return intersec;
248  }
249 
250  //local position
251  intersec = pxdSensorInfo.pointToLocal(gfTrackState.getPos(), m_useAlignment);
252 
253  //try to get the momentum
254  B2DEBUG(1, "Fitted momentum on the plane p = " << gfTrackState.getMom().Mag());
255 
256  // no tolerance currently! Maybe one should be added!
257  double tolerance = 0.0;
258  bool inside = pxdSensorInfo.inside(intersec.X(), intersec.Y(), tolerance, tolerance);
259 
260  // get intersection point in local coordinates with covariance matrix
261  TMatrixDSym covMatrix = gfTrackState.getCov(); // 5D with elements q/p,u',v',u,v in plane system
262 
263  // get ROI by covariance matrix and local intersection point
264  du = std::sqrt(covMatrix(3, 3));
265  dv = std::sqrt(covMatrix(4, 4));
266 
267  if (inside) isgood = true;
268 
269  return intersec;
270 }
271 
272 
273 int
275 {
276  int closest = -1;
277  double mindist = 999999999999; //definitely outside of the sensor
278 
280 
281  //loop the clusters
282  for (int iclus = 0; iclus < m_pxdclusters.getEntries(); iclus++) {
283  //Do not consider as different if only segment differs!
284  //As of this writing segment is never filled for clusters, but just to be sure
285  VxdID clusterID = m_pxdclusters[iclus]->getSensorID();
286  if (avxdid.getLayerNumber() != clusterID.getLayerNumber() ||
287  avxdid.getLadderNumber() != clusterID.getLadderNumber() ||
288  avxdid.getSensorNumber() != clusterID.getSensorNumber()) {
289  continue;
290  }
291  //only cluster on the correct sensor and direction should survive
292 
293  double u = m_pxdclusters[iclus]->getU();
294  double v = m_pxdclusters[iclus]->getV();
295  TVector3 current(u, v, 0);
296 
297  //2D dist sqared
298  double dist = (intersection - current).Mag();
299  if (dist < mindist) {
300  closest = iclus;
301  mindist = dist;
302  }
303  }
304 
305  return closest;
306 
307 }
308 
309 bool PXDDQMEfficiencyNtupleSelftrackModule::isCloseToBorder(int u, int v, int checkDistance)
310 {
311 
312  if (u - checkDistance < 0 || u + checkDistance >= 250 ||
313  v - checkDistance < 0 || v + checkDistance >= 768) {
314  return true;
315  }
316  return false;
317 }
318 
319 bool PXDDQMEfficiencyNtupleSelftrackModule::isDeadPixelClose(int u, int v, int checkDistance, const VxdID& moduleID)
320 {
321 
322  //Iterate over square around the intersection to see if any close pixel is dead
323  for (int u_iter = u - checkDistance; u_iter <= u + checkDistance ; ++u_iter) {
324  for (int v_iter = v - checkDistance; v_iter <= v + checkDistance ; ++v_iter) {
325  if (PXD::PXDPixelMasker::getInstance().pixelDead(moduleID, u_iter, v_iter)
326  || !PXD::PXDPixelMasker::getInstance().pixelOK(moduleID, u_iter, v_iter)) {
327  return true;
328  }
329  }
330  }
331  return false;
332 }
static const ChargedStable pion
charged pion particle
Definition: Const.h:542
Base class for Modules.
Definition: Module.h:72
void initialize() override final
initializes the need store arrays, trees and histograms
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
unsigned int m_minSVDHits
Required hits in SVD strips for tracks.
StoreArray< PXDCluster > m_pxdclusters
store array of pxd clusters
bool isDeadPixelClose(int u, int v, int checkDistance, const VxdID &moduleID)
is a dead pixel close
void terminate() override final
terminate , save tuple to file if needed
void event() override final
main function which fills trees and histograms
std::string m_pxdClustersName
name of the store array of pxd clusters
std::string m_recoTracksName
name of the store array of recotracks
int m_maskedDistance
Distance inside which no dead pixel or module border is allowed.
StoreArray< RecoTrack > m_recoTracks
store array of reco tracks
std::string m_tracksName
name of the store array of tracks
int findClosestCluster(const VxdID &vxdid, TVector3 intersection)
find the closest cluster
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
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
Values of the result of a track fit with a given particle hypothesis.
double getD0() const
Getter for d0.
double getZ0() const
Getter for z0.
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.