Belle II Software development
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#include <framework/geometry/VectorUtil.h>
16
17#include "TDirectory.h"
18#include "TMatrixDSym.h"
19using namespace Belle2;
20
21//-----------------------------------------------------------------
22// Register the Module
23//-----------------------------------------------------------------
24REG_MODULE(PXDDQMEfficiencySelftrack);
25
26//-----------------------------------------------------------------
27// Implementation
28//-----------------------------------------------------------------
29
30PXDDQMEfficiencySelftrackModule::PXDDQMEfficiencySelftrackModule() : HistoModule(), m_vxdGeometry(VXD::GeoCache::getInstance())
31{
32 // Set module properties
33 setDescription("Create basic histograms for PXD efficiency");
34
35 // What exactly is needed for this to be true?
37
38 // Parameter definitions
39 addParam("pxdClustersName", m_pxdClustersName, "name of StoreArray with PXD cluster", std::string(""));
40 addParam("tracksName", m_tracksName, "name of StoreArray with RecoTracks", std::string(""));
41 addParam("ROIsName", m_ROIsName, "name of the list of HLT ROIs, if available in output", 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(4));
45 addParam("binsV", m_v_bins, "histogram bins in v direction", int(6));
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);
76}
77
79{
80 for (auto& h : m_h_track_hits) if (h.second) h.second->Reset();
81 for (auto& h : m_h_matched_cluster) if (h.second) h.second->Reset();
82 for (auto& h : m_h_p) if (h.second) h.second->Reset();
83 for (auto& h : m_h_pt) if (h.second) h.second->Reset();
84 for (auto& h : m_h_su) if (h.second) h.second->Reset();
85 for (auto& h : m_h_sv) if (h.second) h.second->Reset();
86 for (auto& h : m_h_p2) if (h.second) h.second->Reset();
87 for (auto& h : m_h_pt2) if (h.second) h.second->Reset();
88 for (auto& h : m_h_su2) if (h.second) h.second->Reset();
89 for (auto& h : m_h_sv2) if (h.second) h.second->Reset();
90}
91
93{
94 if (!m_pxdclusters.isValid()) {
95 B2INFO("PXDClusters array is missing, no efficiencies");
96 return;
97 }
98 if (!m_tracks.isValid()) {
99 B2INFO("RecoTrack array is missing, no efficiencies");
100 return;
101 }
102 if (!m_ROIs.isValid() && m_requireROIs) {
103 B2INFO("ROI array is missing but required hits in ROIs, aborting");
104 return;
105 }
106
107
108 for (auto& a_track : m_tracks) {
109
110 //If fit failed assume position pointed to is useless anyway
111 if (!a_track.wasFitSuccessful()) continue;
112
113 if (a_track.getNumberOfSVDHits() < m_minSVDHits) continue;
114
115 const genfit::FitStatus* fitstatus = a_track.getTrackFitStatus();
116 if (fitstatus->getPVal() < m_pcut) continue;
117
118 genfit::MeasuredStateOnPlane trackstate;
119 trackstate = a_track.getMeasuredStateOnPlaneFromFirstHit();
120 if (trackstate.getMom().Mag() < m_momCut) continue;
121 if (trackstate.getMom().Pt() < m_pTCut) continue;
122
123 auto ptr = a_track.getRelated<Track>("Tracks");
124
125 if (!ptr) {
126 B2ERROR("expect a track for fitted recotracks");
127 continue;
128 }
129 auto ptr2 = ptr->getTrackFitResultWithClosestMass(Const::pion);
130 if (!ptr2) {
131 B2ERROR("expect a track fit result for mass");
132 continue;
133 }
134
135 // Vertex cut
136 if (ptr2->getZ0() < m_z0minCut || ptr2->getZ0() > m_z0maxCut || fabs(ptr2->getD0()) > m_d0Cut) continue;
137
138 //loop over all PXD sensors to get the intersections
139 std::vector<VxdID> sensors = m_vxdGeometry.getListOfSensors();
140 for (VxdID& aVxdID : sensors) {
142 if (info.getType() != VXD::SensorInfoBase::PXD) continue;
143 //Search for intersections of the track with all PXD layers
144 //Traditional (aka the person before did it like this) method
145 //If there is a way to find out sensors crossed by a track directly, that would most likely be faster
146
147 bool isgood = false;
148 //true = track intersects current sensor
149 double sigu(-9999);
150 double sigv(-9999);
151 ROOT::Math::XYZVector intersec_buff = getTrackInterSec(info, a_track, isgood, sigu, sigv);
152
153 if (!isgood) {
154 continue;//track does not go through this sensor-> nothing to measure anyway
155 } else {
156 if (m_verboseHistos) {
157 if (m_h_p[aVxdID]) m_h_p[aVxdID]->Fill(trackstate.getMom().Mag());
158 if (m_h_pt[aVxdID]) m_h_pt[aVxdID]->Fill(trackstate.getMom().Pt());
159 if (m_h_su[aVxdID]) m_h_su[aVxdID]->Fill(sigu);
160 if (m_h_sv[aVxdID]) m_h_sv[aVxdID]->Fill(sigv);
161 }
162 if (m_uFactor * sigu > m_distcut) continue; // Error ufak*SigmaU > cut
163 if (m_vFactor * sigv > m_distcut) continue; // Error vfak*SigmaV > cut
164
165 double u_fit = intersec_buff.X();
166 double v_fit = intersec_buff.Y();
167
168 int ucell_fit = info.getUCellID(intersec_buff.X());
169 int vcell_fit = info.getVCellID(intersec_buff.Y());
170
171 if (m_cutBorders && isCloseToBorder(ucell_fit, vcell_fit, m_maskedDistance)) {
172 continue;
173 }
174
175 if (m_maskDeadPixels && isDeadPixelClose(ucell_fit, vcell_fit, m_maskedDistance, aVxdID)) {
176 continue;
177 }
178
179 if (m_requireROIs) {
180 //Check if the intersection is inside a ROI
181 //If not, even if measured the cluster was thrown away->Not PXD's fault
182 bool fitInsideROI = false;
183 for (auto& roit : m_ROIs) {
184 if (aVxdID != roit.getSensorID()) {
185 continue; //ROI on other sensor
186 }
187
188 if (ucell_fit < roit.getMaxUid()
189 && ucell_fit > roit.getMinUid()
190 && vcell_fit < roit.getMaxVid()
191 && vcell_fit > roit.getMinVid()) {
192 fitInsideROI = true;
193 }
194 }
195 if (!fitInsideROI) {
196 continue;//Hit wouldn't have been recorded
197 }
198 }
199
200 //This track should be on the sensor
201 m_h_track_hits[aVxdID]->Fill(ucell_fit, vcell_fit);
202
203 //Now check if the sensor measured a hit here
204
205 int bestcluster = findClosestCluster(aVxdID, intersec_buff);
206 if (bestcluster >= 0) {
207 double u_clus = m_pxdclusters[bestcluster]->getU();
208 double v_clus = m_pxdclusters[bestcluster]->getV();
209
210 //is the closest cluster close enough to the track to count as measured?
211 ROOT::Math::XYZVector dist_clus(u_fit - u_clus, v_fit - v_clus, 0);
212 if (dist_clus.R() <= m_distcut) {
213 m_h_matched_cluster[aVxdID]->Fill(ucell_fit, vcell_fit);
214 if (m_verboseHistos) {
215 if (m_h_p2[aVxdID]) m_h_p2[aVxdID]->Fill(trackstate.getMom().Mag());
216 if (m_h_pt2[aVxdID]) m_h_pt2[aVxdID]->Fill(trackstate.getMom().Pt());
217 if (m_h_su2[aVxdID]) m_h_su2[aVxdID]->Fill(sigu);
218 if (m_h_sv2[aVxdID]) m_h_sv2[aVxdID]->Fill(sigv);
219 }
220 }
221 }
222 }
223 }
224 }
225}
226
227
228
229ROOT::Math::XYZVector
231 bool& isgood, double& du, double& dv)
232{
233 //will be set true if the intersect was found
234 isgood = false;
235
236 ROOT::Math::XYZVector intersec(99999999, 9999999, 0); //point outside the sensor
237
238 genfit::MeasuredStateOnPlane gfTrackState = aTrack.getMeasuredStateOnPlaneFromFirstHit();
239
240 //adopted (aka stolen) from tracking/modules/pxdClusterRescue/PXDClusterRescueROIModule
241 try {
242 // get sensor plane
243 ROOT::Math::XYZVector zeroVec(0, 0, 0);
244 ROOT::Math::XYZVector uVec(1, 0, 0);
245 ROOT::Math::XYZVector vVec(0, 1, 0);
246
247 genfit::DetPlane* sensorPlane = new genfit::DetPlane();
248 sensorPlane->setO(XYZToTVector(pxdSensorInfo.pointToGlobal(zeroVec, m_useAlignment)));
249 sensorPlane->setUV(XYZToTVector(pxdSensorInfo.vectorToGlobal(uVec, m_useAlignment)),
250 XYZToTVector(pxdSensorInfo.vectorToGlobal(vVec, m_useAlignment)));
251
252 //boost pointer (will be deleted automatically ?!?!?)
253 genfit::SharedPlanePtr sensorPlaneSptr(sensorPlane);
254
255 // do extrapolation
256 gfTrackState.extrapolateToPlane(sensorPlaneSptr);
257 } catch (genfit::Exception& gfException) {
258 B2WARNING("Fitting failed: " << gfException.getExcString());
259 isgood = false;
260 return intersec;
261 } catch (...) {
262 B2WARNING("Fitting failed: for some reason");
263 isgood = false;
264 return intersec;
265 }
266
267 //local position
268 intersec = pxdSensorInfo.pointToLocal(ROOT::Math::XYZVector(gfTrackState.getPos()), m_useAlignment);
269
270 //try to get the momentum
271 B2DEBUG(1, "Fitted momentum on the plane p = " << gfTrackState.getMom().Mag());
272
273 // no tolerance currently! Maybe one should be added!
274 double tolerance = 0.0;
275 bool inside = pxdSensorInfo.inside(intersec.X(), intersec.Y(), tolerance, tolerance);
276
277 // get intersection point in local coordinates with covariance matrix
278 TMatrixDSym covMatrix = gfTrackState.getCov(); // 5D with elements q/p,u',v',u,v in plane system
279
280 // get ROI by covariance matrix and local intersection point
281 du = std::sqrt(covMatrix(3, 3));
282 dv = std::sqrt(covMatrix(4, 4));
283
284 if (inside) isgood = true;
285
286 return intersec;
287}
288
289
291{
292 // Create a separate histogram directory and cd into it.
293 TDirectory* oldDir = gDirectory;
294 if (m_histogramDirectoryName != "") {
295 oldDir->mkdir(m_histogramDirectoryName.c_str());
296 oldDir->cd(m_histogramDirectoryName.c_str());
297 }
298
299 std::vector<VxdID> sensors = m_vxdGeometry.getListOfSensors();
300 for (VxdID& avxdid : sensors) {
302 if (info.getType() != VXD::SensorInfoBase::PXD) continue;
303 //Only interested in PXD sensors
304
305 TString buff = (std::string)avxdid;
306 buff.ReplaceAll(".", "_");
307
308 int nu = info.getUCells();
309 int nv = info.getVCells();
310
311 //nu + 1,nv + 1 Bin numbers when using one bin per pixel
312 m_h_track_hits[avxdid] = new TH2F("track_hits_" + buff, "tracks through sensor " + buff,
313 m_u_bins, -0.5, nu - 0.5, m_v_bins, -0.5, nv - 0.5);
314 m_h_matched_cluster[avxdid] = new TH2F("matched_cluster_" + buff, "clusters matched to track intersections " + buff,
315 m_u_bins, -0.5, nu - 0.5, m_v_bins, -0.5, nv - 0.5);
316
317 if (m_verboseHistos) {
318 m_h_p[avxdid] = new TH1F("p_" + buff, "p " + buff, 100, 0, 10);
319 m_h_pt[avxdid] = new TH1F("pt_" + buff, "pt " + buff, 100, 0, 10);
320 m_h_su[avxdid] = new TH1F("su_" + buff, "su " + buff, 1000, 0, 1);
321 m_h_sv[avxdid] = new TH1F("sv_" + buff, "sv " + buff, 1000, 0, 1);
322 m_h_p2[avxdid] = new TH1F("p2_" + buff, "p2 " + buff, 100, 0, 10);
323 m_h_pt2[avxdid] = new TH1F("pt2_" + buff, "pt2 " + buff, 100, 0, 10);
324 m_h_su2[avxdid] = new TH1F("su2_" + buff, "su2 " + buff, 1000, 0, 1);
325 m_h_sv2[avxdid] = new TH1F("sv2_" + buff, "sv2 " + buff, 1000, 0, 1);
326 }
327 }
328 // cd back to root directory
329 oldDir->cd();
330}
331
332
333int
334PXDDQMEfficiencySelftrackModule::findClosestCluster(const VxdID& avxdid, ROOT::Math::XYZVector intersection)
335{
336 int closest = -1;
337 double mindist = 999999999999; //definitely outside of the sensor
338
340
341 //loop the clusters
342 for (int iclus = 0; iclus < m_pxdclusters.getEntries(); iclus++) {
343 //Do not consider as different if only segment differs!
344 //As of this writing segment is never filled for clusters, but just to be sure
345 VxdID clusterID = m_pxdclusters[iclus]->getSensorID();
346 if (avxdid.getLayerNumber() != clusterID.getLayerNumber() ||
347 avxdid.getLadderNumber() != clusterID.getLadderNumber() ||
348 avxdid.getSensorNumber() != clusterID.getSensorNumber()) {
349 continue;
350 }
351 //only cluster on the correct sensor and direction should survive
352
353 double u = m_pxdclusters[iclus]->getU();
354 double v = m_pxdclusters[iclus]->getV();
355 ROOT::Math::XYZVector current(u, v, 0);
356
357 //2D dist sqared
358 double dist = (intersection - current).R();
359 if (dist < mindist) {
360 closest = iclus;
361 mindist = dist;
362 }
363 }
364
365 return closest;
366
367}
368
369bool PXDDQMEfficiencySelftrackModule::isCloseToBorder(int u, int v, int checkDistance)
370{
371
372 if (u - checkDistance < 0 || u + checkDistance >= 250 ||
373 v - checkDistance < 0 || v + checkDistance >= 768) {
374 return true;
375 }
376 return false;
377}
378
379bool PXDDQMEfficiencySelftrackModule::isDeadPixelClose(int u, int v, int checkDistance, const VxdID& moduleID)
380{
381
382 //Iterate over square around the intersection to see if any close pixel is dead
383 for (int u_iter = u - checkDistance; u_iter <= u + checkDistance ; ++u_iter) {
384 for (int v_iter = v - checkDistance; v_iter <= v + checkDistance ; ++v_iter) {
385 if (PXD::PXDPixelMasker::getInstance().pixelDead(moduleID, u_iter, v_iter)
386 || !PXD::PXDPixelMasker::getInstance().pixelOK(moduleID, u_iter, v_iter)) {
387 return true;
388 }
389 }
390 }
391 return false;
392}
double R
typedef autogenerated by FFTW
static const ChargedStable pion
charged pion particle
Definition: Const.h:661
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
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
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
int findClosestCluster(const VxdID &vxdid, ROOT::Math::XYZVector intersection)
find the closest cluster
PXDDQMEfficiencySelftrackModule()
Constructor: Sets the description, the properties and the parameters of the module.
void defineHisto() override final
actually defines the trees and histograms
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
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
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: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
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.
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
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
#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.