Belle II Software  release-05-01-25
PXDDQMEfficiencyNtupleModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2013 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Bjoern Spruck *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <pxd/modules/pxdDQM/PXDDQMEfficiencyNtupleModule.h>
12 #include <tracking/dataobjects/ROIid.h>
13 
14 #include <pxd/reconstruction/PXDPixelMasker.h>
15 #include <mdst/dataobjects/Track.h>
16 #include <framework/gearbox/Const.h>
17 
18 #include "TMatrixDSym.h"
19 using namespace Belle2;
20 
21 //-----------------------------------------------------------------
22 // Register the Module
23 //-----------------------------------------------------------------
24 REG_MODULE(PXDDQMEfficiencyNtuple)
25 
26 //-----------------------------------------------------------------
27 // Implementation
28 //-----------------------------------------------------------------
29 
30 PXDDQMEfficiencyNtupleModule::PXDDQMEfficiencyNtupleModule() : Module(), 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("ntupleName", m_ntupleName, "name of ntuple file", std::string("test.root"));
39  addParam("pxdClustersName", m_pxdClustersName, "name of StoreArray with PXD cluster", std::string(""));
40  addParam("recoTracksName", m_recoTracksName, "name of StoreArray with RecoTracks", std::string(""));
41  addParam("tracksName", m_tracksName, "name of StoreArray with Tracks", std::string(""));
42  addParam("ROIsName", m_ROIsName, "name of the list of HLT ROIs, if available in output", std::string(""));
43  addParam("PXDInterceptListName", m_PXDInterceptListName, "name of the list of interceptions", 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(m_ntupleName.data(), "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
83  m_pxdclusters.isOptional(m_pxdClustersName);
84  m_recoTracks.isOptional(m_recoTracksName);
85  m_tracks.isOptional(m_tracksName);
86  m_ROIs.isOptional(m_ROIsName);
88 }
89 
90 
92 {
93  if (!m_pxdclusters.isValid()) {
94  B2INFO("PXDClusters array is missing, no efficiencies");
95  return;
96  }
97  if (!m_recoTracks.isValid()) {
98  B2INFO("RecoTrack array is missing, no efficiencies");
99  return;
100  }
101  if (!m_tracks.isValid()) {
102  B2INFO("Track array is missing, no efficiencies");
103  return;
104  }
105  if (!m_intercepts.isValid()) {
106  B2INFO("Intercept array is missing, no efficiencies");
107  return;
108  }
109 
110  for (auto& track : m_tracks) {
111  RelationVector<RecoTrack> recoTrack = track.getRelationsTo<RecoTrack>(m_recoTracksName);
112  if (!recoTrack.size()) continue;
113 
114  auto a_track = recoTrack[0];
115  //If fit failed assume position pointed to is useless anyway
116  if (!a_track->wasFitSuccessful()) continue;
117 
118  if (a_track->getNumberOfSVDHits() < m_minSVDHits) continue;
119 
120  RelationVector<PXDIntercept> interceptList = a_track->getRelationsTo<PXDIntercept>(m_PXDInterceptListName);
121  if (!interceptList.size()) continue;
122 
123  const genfit::FitStatus* fitstatus = a_track->getTrackFitStatus();
124  if (fitstatus->getPVal() < m_pcut) continue;
125 
126  genfit::MeasuredStateOnPlane trackstate;
127  trackstate = a_track->getMeasuredStateOnPlaneFromFirstHit();
128  if (trackstate.getMom().Mag() < m_momCut) continue;
129  if (trackstate.getMom().Pt() < m_pTCut) continue;
130 
131  const TrackFitResult* ptr2 = track.getTrackFitResultWithClosestMass(Const::pion);
132  if (!ptr2) {
133  B2ERROR("expect a track fit result for mass");
134  continue;
135  }
136 
137  //loop over all PXD sensors to get the intersections
138  std::vector<VxdID> sensors = m_vxdGeometry.getListOfSensors();
139  for (auto intercept : interceptList) {
140  auto aVxdID = intercept.getSensorID();
141  auto& info = m_vxdGeometry.getSensorInfo(aVxdID);
142  //Search for intersections of the track with all PXD layers
143  //Traditional (aka the person before did it like this) method
144  //If there is a way to find out sensors crossed by a track directly, that would most likely be faster
145 
146  double sigu = intercept.getSigmaU();
147  double sigv = intercept.getSigmaV();
148 
149  double u_fit = intercept.getCoorU();
150  double v_fit = intercept.getCoorV();
151 
152  int ucell_fit = info.getUCellID(u_fit); // check wie overflow!!!
153  int vcell_fit = info.getVCellID(v_fit); // Check wie overflow
154 
155  bool closeToBoarder = false;
156  if (isCloseToBorder(ucell_fit, vcell_fit, m_maskedDistance)) {
157  closeToBoarder = true;
158  }
159 
160  bool closeToDead = false;
161  if (isDeadPixelClose(ucell_fit, vcell_fit, m_maskedDistance, aVxdID)) {
162  closeToDead = true;
163  }
164 
165  bool fitInsideROI = false;
166  if (m_ROIs.isValid()) {
167  //Check if the intersection is inside a ROI
168  //If not, even if measured the cluster was thrown away->Not PXD's fault
169  for (auto& roit : m_ROIs) {
170  if (aVxdID != roit.getSensorID()) {
171  continue; //ROI on other sensor
172  }
173 
174  if (ucell_fit < roit.getMaxUid()
175  && ucell_fit > roit.getMinUid()
176  && vcell_fit < roit.getMaxVid()
177  && vcell_fit > roit.getMinVid()) {
178  fitInsideROI = true;
179  }
180  }
181  }
182 
183  //Now check if the sensor measured a hit here
184 
185  int bestcluster = findClosestCluster(aVxdID, TVector3(u_fit, v_fit, 0));
186  double du_clus = 0;
187  double dv_clus = 0;
188  double d_clus = 0;
189  float charge = 0;
190  bool matched = false;
191  if (bestcluster >= 0) {
192  double u_clus = m_pxdclusters[bestcluster]->getU();
193  double v_clus = m_pxdclusters[bestcluster]->getV();
194 
195  //is the closest cluster close enough to the track to count as measured?
196  TVector3 dist_clus(u_fit - u_clus, v_fit - v_clus, 0);
197  du_clus = u_fit - u_clus;
198  dv_clus = v_fit - v_clus;
199  d_clus = dist_clus.Mag();
200  charge = m_pxdclusters[bestcluster]->getCharge();
201  matched = true;
202  }
203  float fill[22] = {float((int)aVxdID), float(u_fit), float(v_fit), float(trackstate.getMom().Mag()), float(trackstate.getMom().Pt()),
204  float(du_clus), float(dv_clus), float(sigu), float(sigv), float(d_clus),
205  float(fitInsideROI), float(closeToBoarder), float(closeToDead), float(matched),
206  float(ptr2->getZ0()), float(ptr2->getD0()), float(a_track->getNumberOfSVDHits()),
207  charge, float(trackstate.getMom().Phi()), float(trackstate.getMom().CosTheta())
208  };
209  m_tuple->Fill(fill);
210  }
211  }
212 }
213 
214 int
215 PXDDQMEfficiencyNtupleModule::findClosestCluster(const VxdID& avxdid, TVector3 intersection)
216 {
217  int closest = -1;
218  double mindist = 999999999999; //definitely outside of the sensor
219 
221 
222  //loop the clusters
223  for (int iclus = 0; iclus < m_pxdclusters.getEntries(); iclus++) {
224  //Do not consider as different if only segment differs!
225  //As of this writing segment is never filled for clusters, but just to be sure
226  VxdID clusterID = m_pxdclusters[iclus]->getSensorID();
227  if (avxdid.getLayerNumber() != clusterID.getLayerNumber() ||
228  avxdid.getLadderNumber() != clusterID.getLadderNumber() ||
229  avxdid.getSensorNumber() != clusterID.getSensorNumber()) {
230  continue;
231  }
232  //only cluster on the correct sensor and direction should survive
233 
234  double u = m_pxdclusters[iclus]->getU();
235  double v = m_pxdclusters[iclus]->getV();
236  TVector3 current(u, v, 0);
237 
238  //2D dist sqared
239  double dist = (intersection - current).Mag();
240  if (dist < mindist) {
241  closest = iclus;
242  mindist = dist;
243  }
244  }
245 
246  return closest;
247 
248 }
249 
250 bool PXDDQMEfficiencyNtupleModule::isCloseToBorder(int u, int v, int checkDistance)
251 {
252 
253  if (u - checkDistance < 0 || u + checkDistance >= 250 ||
254  v - checkDistance < 0 || v + checkDistance >= 768) {
255  return true;
256  }
257  return false;
258 }
259 
260 bool PXDDQMEfficiencyNtupleModule::isDeadPixelClose(int u, int v, int checkDistance, const VxdID& moduleID)
261 {
262 
263  //Iterate over square around the intersection to see if any close pixel is dead
264  for (int u_iter = u - checkDistance; u_iter <= u + checkDistance ; ++u_iter) {
265  for (int v_iter = v - checkDistance; v_iter <= v + checkDistance ; ++v_iter) {
266  if (PXD::PXDPixelMasker::getInstance().pixelDead(moduleID, u_iter, v_iter)
267  || !PXD::PXDPixelMasker::getInstance().pixelOK(moduleID, u_iter, v_iter)) {
268  return true;
269  }
270  }
271  }
272  return false;
273 }
Belle2::RelationVector::size
size_t size() const
Get number of relations.
Definition: RelationVector.h:98
Belle2::PXDDQMEfficiencyNtupleModule::m_file
TFile * m_file
pointer to opened file
Definition: PXDDQMEfficiencyNtupleModule.h:114
Belle2::PXDDQMEfficiencyNtupleModule::isDeadPixelClose
bool isDeadPixelClose(int u, int v, int checkDistance, const VxdID &moduleID)
is a dead pixel close
Definition: PXDDQMEfficiencyNtupleModule.cc:260
Belle2::VxdID
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:43
Belle2::PXDDQMEfficiencyNtupleModule::m_momCut
double m_momCut
Cut on fitted track momentum.
Definition: PXDDQMEfficiencyNtupleModule.h:109
Belle2::PXDDQMEfficiencyNtupleModule::m_pcut
double m_pcut
pValue-Cut for tracks
Definition: PXDDQMEfficiencyNtupleModule.h:108
genfit::MeasuredStateOnPlane
#StateOnPlane with additional covariance matrix.
Definition: MeasuredStateOnPlane.h:39
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::PXDDQMEfficiencyNtupleModule::m_intercepts
StoreArray< PXDIntercept > m_intercepts
store array of PXD Intercepts
Definition: PXDDQMEfficiencyNtupleModule.h:106
Belle2::PXDDQMEfficiencyNtupleModule::m_tracks
StoreArray< Track > m_tracks
store array of tracks
Definition: PXDDQMEfficiencyNtupleModule.h:103
Belle2::PXDDQMEfficiencyNtupleModule::m_pxdClustersName
std::string m_pxdClustersName
name of the store array of pxd clusters
Definition: PXDDQMEfficiencyNtupleModule.h:96
Belle2::PXDDQMEfficiencyNtupleModule::m_tracksName
std::string m_tracksName
name of the store array of tracks
Definition: PXDDQMEfficiencyNtupleModule.h:97
Belle2::VxdID::getLadderNumber
baseType getLadderNumber() const
Get the ladder id.
Definition: VxdID.h:108
Belle2::PXDDQMEfficiencyNtupleModule::event
void event() override final
main function which fills trees and histograms
Definition: PXDDQMEfficiencyNtupleModule.cc:91
Belle2::VXD::SensorInfoBase
Base class to provide Sensor Information for PXD and SVD.
Definition: SensorInfoBase.h:40
Belle2::PXDDQMEfficiencyNtupleModule::m_recoTracksName
std::string m_recoTracksName
name of the store array of recotracks
Definition: PXDDQMEfficiencyNtupleModule.h:98
Belle2::PXDDQMEfficiencyNtupleModule::m_pxdclusters
StoreArray< PXDCluster > m_pxdclusters
store array of pxd clusters
Definition: PXDDQMEfficiencyNtupleModule.h:102
Belle2::PXDDQMEfficiencyNtupleModule::m_minSVDHits
unsigned int m_minSVDHits
Required hits in SVD strips for tracks.
Definition: PXDDQMEfficiencyNtupleModule.h:111
Belle2::TrackFitResult
Values of the result of a track fit with a given particle hypothesis.
Definition: TrackFitResult.h:59
Belle2::VXD::GeoCache::getListOfSensors
const std::vector< VxdID > getListOfSensors() const
Get list of all sensors.
Definition: GeoCache.cc:60
Belle2::TrackFitResult::getZ0
double getZ0() const
Getter for z0.
Definition: TrackFitResult.h:200
Belle2::PXDDQMEfficiencyNtupleModule::initialize
void initialize() override final
initializes the need store arrays, trees and histograms
Definition: PXDDQMEfficiencyNtupleModule.cc:74
Belle2::Module
Base class for Modules.
Definition: Module.h:74
Belle2::RecoTrack
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:78
Belle2::Const::pion
static const ChargedStable pion
charged pion particle
Definition: Const.h:535
genfit::FitStatus::getPVal
virtual double getPVal() const
Get the p value of the fit.
Definition: FitStatus.h:128
Belle2::RelationVector
Class for type safe access to objects that are referred to in relations.
Definition: DataStore.h:38
Belle2::PXDDQMEfficiencyNtupleModule::m_PXDInterceptListName
std::string m_PXDInterceptListName
intercept list name
Definition: PXDDQMEfficiencyNtupleModule.h:100
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::PXDDQMEfficiencyNtupleModule::m_tuple
TNtuple * m_tuple
pointer to opened tuple
Definition: PXDDQMEfficiencyNtupleModule.h:115
Belle2::PXDDQMEfficiencyNtupleModule::m_vxdGeometry
VXD::GeoCache & m_vxdGeometry
the geometry
Definition: PXDDQMEfficiencyNtupleModule.h:93
Belle2::PXDDQMEfficiencyNtupleModule::terminate
void terminate() override final
terminate , save tuple to file if needed
Definition: PXDDQMEfficiencyNtupleModule.cc:53
Belle2::PXDDQMEfficiencyNtupleModule
Creates Ntuples for PXD Efficiency analysis.
Definition: PXDDQMEfficiencyNtupleModule.h:53
Belle2::PXDDQMEfficiencyNtupleModule::isCloseToBorder
bool isCloseToBorder(int u, int v, int checkDistance)
is it close to the border
Definition: PXDDQMEfficiencyNtupleModule.cc:250
Belle2::StoreArray::isValid
bool isValid() const
Check wether the array was registered.
Definition: StoreArray.h:298
Belle2::VxdID::getSensorNumber
baseType getSensorNumber() const
Get the sensor id.
Definition: VxdID.h:110
Belle2::PXDDQMEfficiencyNtupleModule::m_maskedDistance
int m_maskedDistance
Distance inside which no dead pixel or module border is allowed.
Definition: PXDDQMEfficiencyNtupleModule.h:112
Belle2::PXDDQMEfficiencyNtupleModule::m_ROIs
StoreArray< ROIid > m_ROIs
store array of ROIs
Definition: PXDDQMEfficiencyNtupleModule.h:105
Belle2::PXDDQMEfficiencyNtupleModule::m_ntupleName
std::string m_ntupleName
name output file
Definition: PXDDQMEfficiencyNtupleModule.h:95
Belle2::VXD::GeoCache::getSensorInfo
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:68
Belle2::VxdID::getLayerNumber
baseType getLayerNumber() const
Get the layer id.
Definition: VxdID.h:106
Belle2::PXDDQMEfficiencyNtupleModule::m_recoTracks
StoreArray< RecoTrack > m_recoTracks
store array of reco tracks
Definition: PXDDQMEfficiencyNtupleModule.h:104
Belle2::PXDDQMEfficiencyNtupleModule::m_pTCut
double m_pTCut
Cut on fitted track pT.
Definition: PXDDQMEfficiencyNtupleModule.h:110
genfit::FitStatus
Class where important numbers and properties of a fit can be stored.
Definition: FitStatus.h:80
Belle2::PXDIntercept
PXDIntercept stores the U,V coordinates and uncertainties of the intersection of a track with an PXD ...
Definition: PXDIntercept.h:32
Belle2::StoreArray::getEntries
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:226
Belle2::PXDDQMEfficiencyNtupleModule::findClosestCluster
int findClosestCluster(const VxdID &vxdid, TVector3 intersection)
find the closest cluster
Definition: PXDDQMEfficiencyNtupleModule.cc:215
Belle2::PXD::PXDPixelMasker::getInstance
static PXDPixelMasker & getInstance()
Main (and only) way to access the PXDPixelMasker.
Definition: PXDPixelMasker.cc:37
Belle2::TrackFitResult::getD0
double getD0() const
Getter for d0.
Definition: TrackFitResult.h:178
Belle2::PXDDQMEfficiencyNtupleModule::m_ROIsName
std::string m_ROIsName
name of the store array of ROIs
Definition: PXDDQMEfficiencyNtupleModule.h:99