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