Belle II Software development
PXDDQMEfficiencyNtupleModule.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/PXDDQMEfficiencyNtupleModule.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"
17using namespace Belle2;
18
19//-----------------------------------------------------------------
20// Register the Module
21//-----------------------------------------------------------------
22REG_MODULE(PXDDQMEfficiencyNtuple);
23
24//-----------------------------------------------------------------
25// Implementation
26//-----------------------------------------------------------------
27
28PXDDQMEfficiencyNtupleModule::PXDDQMEfficiencyNtupleModule() : Module(), m_vxdGeometry(VXD::GeoCache::getInstance())
29{
30 // Set module properties
31 setDescription("Create basic histograms for PXD efficiency");
32
33 // setPropertyFlags(c_ParallelProcessingCertified);// for ntuple not certified!!!
34
35 // Parameter definitions
36 addParam("ntupleName", m_ntupleName, "name of ntuple file", std::string("test.root"));
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("PXDInterceptListName", m_PXDInterceptListName, "name of the list of interceptions", 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(m_ntupleName.data(), "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);
86}
87
88
90{
91 if (!m_pxdclusters.isValid()) {
92 B2INFO("PXDClusters array is missing, no efficiencies");
93 return;
94 }
95 if (!m_recoTracks.isValid()) {
96 B2INFO("RecoTrack array is missing, no efficiencies");
97 return;
98 }
99 if (!m_tracks.isValid()) {
100 B2INFO("Track array is missing, no efficiencies");
101 return;
102 }
103 if (!m_intercepts.isValid()) {
104 B2INFO("Intercept array is missing, no efficiencies");
105 return;
106 }
107
108 for (auto& track : m_tracks) {
109 RelationVector<RecoTrack> recoTrack = track.getRelationsTo<RecoTrack>(m_recoTracksName);
110 if (!recoTrack.size()) continue;
111
112 auto a_track = recoTrack[0];
113 //If fit failed assume position pointed to is useless anyway
114 if (!a_track->wasFitSuccessful()) continue;
115
116 if (a_track->getNumberOfSVDHits() < m_minSVDHits) continue;
117
118 RelationVector<PXDIntercept> interceptList = a_track->getRelationsTo<PXDIntercept>(m_PXDInterceptListName);
119 if (!interceptList.size()) continue;
120
121 const genfit::FitStatus* fitstatus = a_track->getTrackFitStatus();
122 if (fitstatus->getPVal() < m_pcut) continue;
123
124 genfit::MeasuredStateOnPlane trackstate;
125 trackstate = a_track->getMeasuredStateOnPlaneFromFirstHit();
126 if (trackstate.getMom().Mag() < m_momCut) continue;
127 if (trackstate.getMom().Pt() < m_pTCut) continue;
128
129 const TrackFitResult* ptr2 = track.getTrackFitResultWithClosestMass(Const::pion);
130 if (!ptr2) {
131 B2ERROR("expect a track fit result for mass");
132 continue;
133 }
134
135 //loop over all PXD sensors to get the intersections
136 std::vector<VxdID> sensors = m_vxdGeometry.getListOfSensors();
137 for (auto intercept : interceptList) {
138 auto aVxdID = intercept.getSensorID();
139 auto& info = m_vxdGeometry.getSensorInfo(aVxdID);
140 //Search for intersections of the track with all PXD layers
141 //Traditional (aka the person before did it like this) method
142 //If there is a way to find out sensors crossed by a track directly, that would most likely be faster
143
144 double sigu = intercept.getSigmaU();
145 double sigv = intercept.getSigmaV();
146
147 double u_fit = intercept.getCoorU();
148 double v_fit = intercept.getCoorV();
149
150 int ucell_fit = info.getUCellID(u_fit); // check wie overflow!!!
151 int vcell_fit = info.getVCellID(v_fit); // Check wie overflow
152
153 bool closeToBoarder = false;
154 if (isCloseToBorder(ucell_fit, vcell_fit, m_maskedDistance)) {
155 closeToBoarder = true;
156 }
157
158 bool closeToDead = false;
159 if (isDeadPixelClose(ucell_fit, vcell_fit, m_maskedDistance, aVxdID)) {
160 closeToDead = true;
161 }
162
163 bool fitInsideROI = false;
164 if (m_ROIs.isValid()) {
165 //Check if the intersection is inside a ROI
166 //If not, even if measured the cluster was thrown away->Not PXD's fault
167 for (auto& roit : m_ROIs) {
168 if (aVxdID != (roit.getSensorID()).getID()) {
169 continue; //ROI on other sensor
170 }
171
172 if (ucell_fit < roit.getMaxUid()
173 && ucell_fit > roit.getMinUid()
174 && vcell_fit < roit.getMaxVid()
175 && vcell_fit > roit.getMinVid()) {
176 fitInsideROI = true;
177 }
178 }
179 }
180
181 //Now check if the sensor measured a hit here
182
183 int bestcluster = findClosestCluster(aVxdID, ROOT::Math::XYZVector(u_fit, v_fit, 0));
184 double du_clus = 0;
185 double dv_clus = 0;
186 double d_clus = 0;
187 float charge = 0;
188 bool matched = false;
189 if (bestcluster >= 0) {
190 double u_clus = m_pxdclusters[bestcluster]->getU();
191 double v_clus = m_pxdclusters[bestcluster]->getV();
192
193 //is the closest cluster close enough to the track to count as measured?
194 ROOT::Math::XYZVector dist_clus(u_fit - u_clus, v_fit - v_clus, 0);
195 du_clus = u_fit - u_clus;
196 dv_clus = v_fit - v_clus;
197 d_clus = dist_clus.R();
198 charge = m_pxdclusters[bestcluster]->getCharge();
199 matched = true;
200 }
201 float fill[22] = {float((int)aVxdID), float(u_fit), float(v_fit), float(trackstate.getMom().Mag()), float(trackstate.getMom().Pt()),
202 float(du_clus), float(dv_clus), float(sigu), float(sigv), float(d_clus),
203 float(fitInsideROI), float(closeToBoarder), float(closeToDead), float(matched),
204 float(ptr2->getZ0()), float(ptr2->getD0()), float(a_track->getNumberOfSVDHits()),
205 charge, float(trackstate.getMom().Phi()), float(trackstate.getMom().CosTheta())
206 };
207 m_tuple->Fill(fill);
208 }
209 }
210}
211
212int
213PXDDQMEfficiencyNtupleModule::findClosestCluster(const VxdID& avxdid, ROOT::Math::XYZVector intersection)
214{
215 int closest = -1;
216 double mindist = 999999999999; //definitely outside of the sensor
217
219
220 //loop the clusters
221 for (int iclus = 0; iclus < m_pxdclusters.getEntries(); iclus++) {
222 //Do not consider as different if only segment differs!
223 //As of this writing segment is never filled for clusters, but just to be sure
224 VxdID clusterID = m_pxdclusters[iclus]->getSensorID();
225 if (avxdid.getLayerNumber() != clusterID.getLayerNumber() ||
226 avxdid.getLadderNumber() != clusterID.getLadderNumber() ||
227 avxdid.getSensorNumber() != clusterID.getSensorNumber()) {
228 continue;
229 }
230 //only cluster on the correct sensor and direction should survive
231
232 double u = m_pxdclusters[iclus]->getU();
233 double v = m_pxdclusters[iclus]->getV();
234 ROOT::Math::XYZVector current(u, v, 0);
235
236 //2D dist sqared
237 double dist = (intersection - current).R();
238 if (dist < mindist) {
239 closest = iclus;
240 mindist = dist;
241 }
242 }
243
244 return closest;
245
246}
247
248bool PXDDQMEfficiencyNtupleModule::isCloseToBorder(int u, int v, int checkDistance)
249{
250
251 if (u - checkDistance < 0 || u + checkDistance >= 250 ||
252 v - checkDistance < 0 || v + checkDistance >= 768) {
253 return true;
254 }
255 return false;
256}
257
258bool PXDDQMEfficiencyNtupleModule::isDeadPixelClose(int u, int v, int checkDistance, const VxdID& moduleID)
259{
260
261 //Iterate over square around the intersection to see if any close pixel is dead
262 for (int u_iter = u - checkDistance; u_iter <= u + checkDistance ; ++u_iter) {
263 for (int v_iter = v - checkDistance; v_iter <= v + checkDistance ; ++v_iter) {
264 if (PXD::PXDPixelMasker::getInstance().pixelDead(moduleID, u_iter, v_iter)
265 || !PXD::PXDPixelMasker::getInstance().pixelOK(moduleID, u_iter, v_iter)) {
266 return true;
267 }
268 }
269 }
270 return false;
271}
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
std::string m_PXDInterceptListName
intercept list name
double m_momCut
Cut on fitted track momentum.
int findClosestCluster(const VxdID &vxdid, ROOT::Math::XYZVector intersection)
find the closest cluster
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
StoreArray< Track > m_tracks
store array of tracks
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< PXDIntercept > m_intercepts
store array of PXD Intercepts
std::string m_ROIsName
name of the store array of ROIs
StoreArray< RecoTrack > m_recoTracks
store array of reco tracks
bool m_useAlignment
if true alignment will be used!
PXDDQMEfficiencyNtupleModule()
Constructor: Sets the description, the properties and the parameters of the module.
std::string m_tracksName
name of the store array of tracks
StoreArray< ROIid > m_ROIs
store array of ROIs
PXDIntercept stores the U,V coordinates and uncertainties of the intersection of a track with an PXD ...
Definition: PXDIntercept.h:22
static PXDPixelMasker & getInstance()
Main (and only) way to access the PXDPixelMasker.
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:79
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.
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
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.