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