CMS 3D CMS Logo

SiPixelHitEfficiencySource.cc
Go to the documentation of this file.
1 // Package: SiPixelMonitorTrack
2 // Class: SiPixelHitEfficiencySource
3 //
4 // class SiPixelHitEfficiencyModule SiPixelHitEfficiencyModule.cc
5 // DQM/SiPixelMonitorTrack/src/SiPixelHitEfficiencyModule.cc
6 //
7 // Description: SiPixel hit efficiency data quality monitoring modules
8 // Implementation: prototype -> improved -> never final - end of the 1st step
9 //
10 // Original Authors: Romain Rougny & Luca Mucibello
11 // Created: Mar Nov 10 13:29:00 CET 2009
12 
13 #include <cmath>
14 #include <iostream>
15 #include <map>
16 #include <string>
17 #include <utility>
18 #include <vector>
19 
22 
33 
45 
48 
49 using namespace std;
50 using namespace edm;
51 
53  : pSet_(pSet),
54  modOn(pSet.getUntrackedParameter<bool>("modOn", true)),
55  ladOn(pSet.getUntrackedParameter<bool>("ladOn", false)),
56  layOn(pSet.getUntrackedParameter<bool>("layOn", false)),
57  phiOn(pSet.getUntrackedParameter<bool>("phiOn", false)),
58  ringOn(pSet.getUntrackedParameter<bool>("ringOn", false)),
59  bladeOn(pSet.getUntrackedParameter<bool>("bladeOn", false)),
60  diskOn(pSet.getUntrackedParameter<bool>("diskOn", false)),
61  isUpgrade(pSet.getUntrackedParameter<bool>("isUpgrade", false))
62 // updateEfficiencies(
63 // pSet.getUntrackedParameter<bool>("updateEfficiencies",false) )
64 {
65  pSet_ = pSet;
66  debug_ = pSet_.getUntrackedParameter<bool>("debug", false);
67  applyEdgeCut_ = pSet_.getUntrackedParameter<bool>("applyEdgeCut");
68  nSigma_EdgeCut_ = pSet_.getUntrackedParameter<double>("nSigma_EdgeCut");
69  vtxsrc_ = pSet_.getUntrackedParameter<std::string>("vtxsrc", "offlinePrimaryVertices");
70  vertexCollectionToken_ = consumes<reco::VertexCollection>(vtxsrc_);
71  tracksrc_ = consumes<TrajTrackAssociationCollection>(pSet_.getParameter<edm::InputTag>("trajectoryInput"));
72  clusterCollectionToken_ = consumes<edmNew::DetSetVector<SiPixelCluster>>(std::string("siPixelClusters"));
73 
74  measurementTrackerEventToken_ = consumes<MeasurementTrackerEvent>(std::string("MeasurementTrackerEvent"));
75 
76  trackerTopoToken_ = esConsumes<TrackerTopology, TrackerTopologyRcd>();
77  trackerGeomToken_ = esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>();
78  measurementTrackerToken_ = esConsumes<MeasurementTracker, CkfComponentsRecord>();
80  esConsumes<Chi2MeasurementEstimatorBase, TrackingComponentsRecord>(edm::ESInputTag("", "Chi2"));
81  propagatorToken_ = esConsumes<Propagator, TrackingComponentsRecord>(edm::ESInputTag("", "PropagatorWithMaterial"));
83  esConsumes<PixelClusterParameterEstimator, TkPixelCPERecord>(edm::ESInputTag("", "PixelCPEGeneric"));
84  trackerTopoTokenBeginRun_ = esConsumes<TrackerTopology, TrackerTopologyRcd, edm::Transition::BeginRun>();
85  trackerGeomTokenBeginRun_ = esConsumes<TrackerGeometry, TrackerDigiGeometryRecord, edm::Transition::BeginRun>();
86 
87  firstRun = true;
88 
89  LogInfo("PixelDQM") << "SiPixelHitEfficiencySource constructor" << endl;
90  LogInfo("PixelDQM") << "Mod/Lad/Lay/Phi " << modOn << "/" << ladOn << "/" << layOn << "/" << phiOn << std::endl;
91  LogInfo("PixelDQM") << "Blade/Disk/Ring" << bladeOn << "/" << diskOn << "/" << ringOn << std::endl;
92 }
93 
95  LogInfo("PixelDQM") << "SiPixelHitEfficiencySource destructor" << endl;
96 
97  std::map<uint32_t, SiPixelHitEfficiencyModule *>::iterator struct_iter;
98  for (struct_iter = theSiPixelStructure.begin(); struct_iter != theSiPixelStructure.end(); struct_iter++) {
99  delete struct_iter->second;
100  struct_iter->second = nullptr;
101  }
102 }
103 
104 void SiPixelHitEfficiencySource::fillClusterProbability(int layer, int disk, bool plus, double probability) {
105  // barrel
106  if (layer != 0) {
107  if (layer == 1) {
108  if (plus)
109  meClusterProbabilityL1_Plus_->Fill(probability);
110  else
111  meClusterProbabilityL1_Minus_->Fill(probability);
112  }
113 
114  else if (layer == 2) {
115  if (plus)
116  meClusterProbabilityL2_Plus_->Fill(probability);
117  else
118  meClusterProbabilityL2_Minus_->Fill(probability);
119  }
120 
121  else if (layer == 3) {
122  if (plus)
123  meClusterProbabilityL3_Plus_->Fill(probability);
124  else
125  meClusterProbabilityL3_Minus_->Fill(probability);
126  }
127  }
128  // Endcap
129  if (disk != 0) {
130  if (disk == 1) {
131  if (plus)
132  meClusterProbabilityD1_Plus_->Fill(probability);
133  else
134  meClusterProbabilityD1_Minus_->Fill(probability);
135  }
136  if (disk == 2) {
137  if (plus)
138  meClusterProbabilityD2_Plus_->Fill(probability);
139  else
140  meClusterProbabilityD2_Minus_->Fill(probability);
141  }
142  }
143 }
144 
146  LogInfo("PixelDQM") << "SiPixelHitEfficiencySource beginRun()" << endl;
147 
148  if (firstRun) {
149  // retrieve TrackerGeometry for pixel dets
150 
151  nvalid = 0;
152  nmissing = 0;
153 
154  firstRun = false;
155  }
156 
158  if (debug_)
159  LogVerbatim("PixelDQM") << "TrackerGeometry " << &(*TG) << " size is " << TG->dets().size() << endl;
160 
161  // build theSiPixelStructure with the pixel barrel and endcap dets from
162  // TrackerGeometry
163  for (TrackerGeometry::DetContainer::const_iterator pxb = TG->detsPXB().begin(); pxb != TG->detsPXB().end(); pxb++) {
164  if (dynamic_cast<PixelGeomDetUnit const *>((*pxb)) != nullptr) {
165  SiPixelHitEfficiencyModule *module = new SiPixelHitEfficiencyModule((*pxb)->geographicalId().rawId());
166  theSiPixelStructure.insert(
167  pair<uint32_t, SiPixelHitEfficiencyModule *>((*pxb)->geographicalId().rawId(), module));
168  }
169  }
170  for (TrackerGeometry::DetContainer::const_iterator pxf = TG->detsPXF().begin(); pxf != TG->detsPXF().end(); pxf++) {
171  if (dynamic_cast<PixelGeomDetUnit const *>((*pxf)) != nullptr) {
172  SiPixelHitEfficiencyModule *module = new SiPixelHitEfficiencyModule((*pxf)->geographicalId().rawId());
173  theSiPixelStructure.insert(
174  pair<uint32_t, SiPixelHitEfficiencyModule *>((*pxf)->geographicalId().rawId(), module));
175  }
176  }
177  LogInfo("PixelDQM") << "SiPixelStructure size is " << theSiPixelStructure.size() << endl;
178 }
179 
181  edm::Run const &iRun,
182  edm::EventSetup const &iSetup) {
183  // book residual histograms in theSiPixelFolder - one (x,y) pair of histograms
184  // per det
185  SiPixelFolderOrganizer theSiPixelFolder(false);
186 
188  const TrackerTopology *pTT = tTopoHandle.product();
189 
190  for (std::map<uint32_t, SiPixelHitEfficiencyModule *>::iterator pxd = theSiPixelStructure.begin();
191  pxd != theSiPixelStructure.end();
192  pxd++) {
193  if (modOn) {
194  if (theSiPixelFolder.setModuleFolder(iBooker, (*pxd).first, 0, isUpgrade))
195  (*pxd).second->book(pSet_, pTT, iBooker, 0, isUpgrade);
196  else
197  throw cms::Exception("LogicError") << "SiPixelHitEfficiencySource Folder Creation Failed! ";
198  }
199  if (ladOn) {
200  if (theSiPixelFolder.setModuleFolder(iBooker, (*pxd).first, 1, isUpgrade))
201  (*pxd).second->book(pSet_, pTT, iBooker, 1, isUpgrade);
202  else
203  throw cms::Exception("LogicError") << "SiPixelHitEfficiencySource ladder Folder Creation Failed! ";
204  }
205  if (layOn) {
206  if (theSiPixelFolder.setModuleFolder(iBooker, (*pxd).first, 2, isUpgrade))
207  (*pxd).second->book(pSet_, pTT, iBooker, 2, isUpgrade);
208  else
209  throw cms::Exception("LogicError") << "SiPixelHitEfficiencySource layer Folder Creation Failed! ";
210  }
211  if (phiOn) {
212  if (theSiPixelFolder.setModuleFolder(iBooker, (*pxd).first, 3, isUpgrade))
213  (*pxd).second->book(pSet_, pTT, iBooker, 3, isUpgrade);
214  else
215  throw cms::Exception("LogicError") << "SiPixelHitEfficiencySource phi Folder Creation Failed! ";
216  }
217  if (bladeOn) {
218  if (theSiPixelFolder.setModuleFolder(iBooker, (*pxd).first, 4, isUpgrade))
219  (*pxd).second->book(pSet_, pTT, iBooker, 4, isUpgrade);
220  else
221  throw cms::Exception("LogicError") << "SiPixelHitEfficiencySource Blade Folder Creation Failed! ";
222  }
223  if (diskOn) {
224  if (theSiPixelFolder.setModuleFolder(iBooker, (*pxd).first, 5, isUpgrade))
225  (*pxd).second->book(pSet_, pTT, iBooker, 5, isUpgrade);
226  else
227  throw cms::Exception("LogicError") << "SiPixelHitEfficiencySource Disk Folder Creation Failed! ";
228  }
229  if (ringOn) {
230  if (theSiPixelFolder.setModuleFolder(iBooker, (*pxd).first, 6, isUpgrade))
231  (*pxd).second->book(pSet_, pTT, iBooker, 6, isUpgrade);
232  else
233  throw cms::Exception("LogicError") << "SiPixelHitEfficiencySource Ring Folder Creation Failed! ";
234  }
235  }
236 
237  // book cluster probability histos for Barrel and Endcap
238  iBooker.setCurrentFolder("Pixel/Barrel");
239 
241  iBooker.book1D("ClusterProbabilityXY_Layer1_Plus", "ClusterProbabilityXY_Layer1_Plus", 500, -14, 0.1);
242  meClusterProbabilityL1_Plus_->setAxisTitle("Log(ClusterProbability)", 1);
243 
245  iBooker.book1D("ClusterProbabilityXY_Layer1_Minus", "ClusterProbabilityXY_Layer1_Minus", 500, -14, 0.1);
246  meClusterProbabilityL1_Minus_->setAxisTitle("Log(ClusterProbability)", 1);
247 
249  iBooker.book1D("ClusterProbabilityXY_Layer2_Plus", "ClusterProbabilityXY_Layer2_Plus", 500, -14, 0.1);
250  meClusterProbabilityL2_Plus_->setAxisTitle("Log(ClusterProbability)", 1);
251 
253  iBooker.book1D("ClusterProbabilityXY_Layer2_Minus", "ClusterProbabilityXY_Layer2_Minus", 500, -14, 0.1);
254  meClusterProbabilityL2_Minus_->setAxisTitle("Log(ClusterProbability)", 1);
255 
257  iBooker.book1D("ClusterProbabilityXY_Layer3_Plus", "ClusterProbabilityXY_Layer3_Plus", 500, -14, 0.1);
258  meClusterProbabilityL3_Plus_->setAxisTitle("Log(ClusterProbability)", 1);
259 
261  iBooker.book1D("ClusterProbabilityXY_Layer3_Minus", "ClusterProbabilityXY_Layer3_Minus", 500, -14, 0.1);
262  meClusterProbabilityL3_Minus_->setAxisTitle("Log(ClusterProbability)", 1);
263 
264  iBooker.setCurrentFolder("Pixel/Endcap");
265 
267  iBooker.book1D("ClusterProbabilityXY_Disk1_Plus", "ClusterProbabilityXY_Disk1_Plus", 500, -14, 0.1);
268  meClusterProbabilityD1_Plus_->setAxisTitle("Log(ClusterProbability)", 1);
269 
271  iBooker.book1D("ClusterProbabilityXY_Disk1_Minus", "ClusterProbabilityXY_Disk1_Minus", 500, -14, 0.1);
272  meClusterProbabilityD1_Minus_->setAxisTitle("Log(ClusterProbability)", 1);
273 
275  iBooker.book1D("ClusterProbabilityXY_Disk2_Plus", "ClusterProbabilityXY_Disk2_Plus", 500, -14, 0.1);
276  meClusterProbabilityD2_Plus_->setAxisTitle("Log(ClusterProbability)", 1);
277 
279  iBooker.book1D("ClusterProbabilityXY_Disk2_Minus", "ClusterProbabilityXY_Disk2_Minus", 500, -14, 0.1);
280  meClusterProbabilityD2_Minus_->setAxisTitle("Log(ClusterProbability)", 1);
281 }
282 
285  const TrackerTopology *pTT = tTopoHandle.product();
286 
287  edm::Handle<reco::VertexCollection> vertexCollectionHandle;
288  iEvent.getByToken(vertexCollectionToken_, vertexCollectionHandle);
289  if (!vertexCollectionHandle.isValid())
290  return;
291  nvtx_ = 0;
292  vtxntrk_ = -9999;
293  vtxD0_ = -9999.;
294  vtxX_ = -9999.;
295  vtxY_ = -9999.;
296  vtxZ_ = -9999.;
297  vtxndof_ = -9999.;
298  vtxchi2_ = -9999.;
299  const reco::VertexCollection &vertices = *vertexCollectionHandle.product();
300  reco::VertexCollection::const_iterator bestVtx = vertices.end();
301  for (reco::VertexCollection::const_iterator it = vertices.begin(); it != vertices.end(); ++it) {
302  if (!it->isValid())
303  continue;
304  if (vtxntrk_ == -9999 || vtxntrk_ < int(it->tracksSize()) ||
305  (vtxntrk_ == int(it->tracksSize()) && fabs(vtxZ_) > fabs(it->z()))) {
306  vtxntrk_ = it->tracksSize();
307  vtxD0_ = it->position().rho();
308  vtxX_ = it->x();
309  vtxY_ = it->y();
310  vtxZ_ = it->z();
311  vtxndof_ = it->ndof();
312  vtxchi2_ = it->chi2();
313  bestVtx = it;
314  }
315  if (fabs(it->z()) <= 20. && fabs(it->position().rho()) <= 2. && it->ndof() > 4)
316  nvtx_++;
317  }
318  if (nvtx_ < 1)
319  return;
320 
321  // get the map
323  iEvent.getByToken(tracksrc_, match);
324 
325  if (debug_) {
327  labelsForToken(tracksrc_, labels);
328  std::cout << "Trajectories from " << labels.module << std::endl;
329  }
330 
331  if (!match.isValid())
332  return;
333 
334  const TrajTrackAssociationCollection ttac = *(match.product());
335 
336  if (debug_) {
337  std::cout << "+++ NEW EVENT +++" << std::endl;
338  std::cout << "Map entries \t : " << ttac.size() << std::endl;
339  }
340 
341  std::set<SiPixelCluster> clusterSet;
342  // TrajectoryStateCombiner tsoscomb;
343  // define variables for extrapolation
344  int extrapolateFrom_ = 2;
345  int extrapolateTo_ = 1;
346  float maxlxmatch_ = 0.2;
347  float maxlymatch_ = 0.2;
348  bool keepOriginalMissingHit_ = true;
349  edm::ESHandle<MeasurementTracker> measurementTrackerHandle = iSetup.getHandle(measurementTrackerToken_);
350 
352  edm::Handle<MeasurementTrackerEvent> measurementTrackerEventHandle;
353  iEvent.getByToken(measurementTrackerEventToken_, measurementTrackerEventHandle);
355  Propagator *thePropagator = prop.product()->clone();
356  // determines direction of the propagator => inward
357  if (extrapolateFrom_ >= extrapolateTo_) {
359  }
360  TrajectoryStateCombiner trajStateComb;
361  bool debug_ = false;
362 
363  // Loop over track collection
364  for (TrajTrackAssociationCollection::const_iterator it = ttac.begin(); it != ttac.end(); ++it) {
365  // define vector to save extrapolated tracks
366  std::vector<TrajectoryMeasurement> expTrajMeasurements;
367 
368  const edm::Ref<std::vector<Trajectory>> traj_iterator = it->key;
369  // Trajectory Map, extract Trajectory for this track
370  reco::TrackRef trackref = it->val;
371  // tracks++;
372 
373  bool isBpixtrack = false, isFpixtrack = false;
374  int nStripHits = 0;
375  int L1hits = 0;
376  int L2hits = 0;
377  int L3hits = 0;
378  int D1hits = 0;
379  int D2hits = 0;
380  std::vector<TrajectoryMeasurement> tmeasColl = traj_iterator->measurements();
381  std::vector<TrajectoryMeasurement>::const_iterator tmeasIt;
382  // loop on measurements to find out what kind of hits there are
383  for (tmeasIt = tmeasColl.begin(); tmeasIt != tmeasColl.end(); tmeasIt++) {
384  // if(! tmeasIt->updatedState().isValid()) continue; NOT NECESSARY (I
385  // HOPE)
386  TransientTrackingRecHit::ConstRecHitPointer testhit = tmeasIt->recHit();
387  if (testhit->geographicalId().det() != DetId::Tracker)
388  continue;
389  uint testSubDetID = (testhit->geographicalId().subdetId());
390  const DetId &hit_detId = testhit->geographicalId();
391  int hit_layer = 0;
392  int hit_ladder = 0;
393  int hit_mod = 0;
394  int hit_disk = 0;
395 
396  if (testSubDetID == PixelSubdetector::PixelBarrel) {
397  isBpixtrack = true;
398  hit_layer = PixelBarrelName(hit_detId, pTT, isUpgrade).layerName();
399 
400  hit_ladder = PXBDetId(hit_detId).ladder();
401  hit_mod = PXBDetId(hit_detId).module();
402 
403  if (hit_layer == 1)
404  L1hits++;
405  if (hit_layer == 2)
406  L2hits++;
407  if (hit_layer == 3)
408  L3hits++;
409  }
410  if (testSubDetID == PixelSubdetector::PixelEndcap) {
411  isFpixtrack = true;
412  hit_disk = PixelEndcapName(hit_detId, pTT, isUpgrade).diskName();
413 
414  if (hit_disk == 1)
415  D1hits++;
416  if (hit_disk == 2)
417  D2hits++;
418  }
419  if (testSubDetID == StripSubdetector::TIB)
420  nStripHits++;
421  if (testSubDetID == StripSubdetector::TOB)
422  nStripHits++;
423  if (testSubDetID == StripSubdetector::TID)
424  nStripHits++;
425  if (testSubDetID == StripSubdetector::TEC)
426  nStripHits++;
427  // check if last valid hit is in Layer 2 or Disk 1
428 
429  bool lastValidL2 = false;
430  if ((testSubDetID == PixelSubdetector::PixelBarrel && hit_layer == extrapolateFrom_) ||
431  (testSubDetID == PixelSubdetector::PixelEndcap && hit_disk == 1)) {
432  if (testhit->isValid()) {
433  if (tmeasIt == tmeasColl.end() - 1) {
434  lastValidL2 = true;
435  } else {
436  tmeasIt++;
437  TransientTrackingRecHit::ConstRecHitPointer nextRecHit = tmeasIt->recHit();
438  uint nextSubDetID = (nextRecHit->geographicalId().subdetId());
439  int nextlayer = PixelBarrelName(nextRecHit->geographicalId()).layerName();
440  if (nextSubDetID == PixelSubdetector::PixelBarrel && nextlayer == extrapolateTo_) {
441  lastValidL2 = true; //&& !nextRecHit->isValid()) lastValidL2=true;
442  }
443  tmeasIt--;
444  }
445  }
446  } // end check last valid layer
447  if (lastValidL2) {
448  std::vector<const BarrelDetLayer *> pxbLayers =
449  measurementTrackerHandle->geometricSearchTracker()->pixelBarrelLayers();
450  const DetLayer *pxb1 = pxbLayers[extrapolateTo_ - 1];
451  const MeasurementEstimator *estimator = est.product();
452  const LayerMeasurements *theLayerMeasurements =
453  new LayerMeasurements(*measurementTrackerHandle, *measurementTrackerEventHandle);
454  const TrajectoryStateOnSurface tsosPXB2 = tmeasIt->updatedState();
455  expTrajMeasurements = theLayerMeasurements->measurements(*pxb1, tsosPXB2, *thePropagator, *estimator);
456  delete theLayerMeasurements;
457  if (!expTrajMeasurements.empty()) {
458  for (uint p = 0; p < expTrajMeasurements.size(); p++) {
459  TrajectoryMeasurement pxb1TM(expTrajMeasurements[p]);
460  const auto &pxb1Hit = pxb1TM.recHit();
461  // remove hits with rawID == 0
462  if (pxb1Hit->geographicalId().rawId() == 0) {
463  expTrajMeasurements.erase(expTrajMeasurements.begin() + p);
464  continue;
465  }
466  }
467  }
468  //
469  }
470  // check if extrapolated hit to layer 1 one matches the original hit
471  TrajectoryStateOnSurface chkPredTrajState =
472  trajStateComb(tmeasIt->forwardPredictedState(), tmeasIt->backwardPredictedState());
473  if (!chkPredTrajState.isValid())
474  continue;
475  float chkx = chkPredTrajState.globalPosition().x();
476  float chky = chkPredTrajState.globalPosition().y();
477  float chkz = chkPredTrajState.globalPosition().z();
478  LocalPoint chklp = chkPredTrajState.localPosition();
479  if (testSubDetID == PixelSubdetector::PixelBarrel && hit_layer == extrapolateTo_) {
480  // Here we will drop the extrapolated hits if there is a hit and use
481  // that hit
482  vector<int> imatches;
483  size_t imatch = 0;
484  float glmatch = 9999.;
485  for (size_t iexp = 0; iexp < expTrajMeasurements.size(); iexp++) {
486  const DetId &exphit_detId = expTrajMeasurements[iexp].recHit()->geographicalId();
487  int exphit_ladder = PXBDetId(exphit_detId).ladder();
488  int exphit_mod = PXBDetId(exphit_detId).module();
489  int dladder = abs(exphit_ladder - hit_ladder);
490  if (dladder > 10)
491  dladder = 20 - dladder;
492  int dmodule = abs(exphit_mod - hit_mod);
493  if (dladder != 0 || dmodule != 0) {
494  continue;
495  }
496 
497  TrajectoryStateOnSurface predTrajState = expTrajMeasurements[iexp].updatedState();
498  float x = predTrajState.globalPosition().x();
499  float y = predTrajState.globalPosition().y();
500  float z = predTrajState.globalPosition().z();
501  float dxyz = sqrt((chkx - x) * (chkx - x) + (chky - y) * (chky - y) + (chkz - z) * (chkz - z));
502 
503  if (dxyz <= glmatch) {
504  glmatch = dxyz;
505  imatch = iexp;
506  imatches.push_back(int(imatch));
507  }
508 
509  } // found the propagated traj best matching the hit in data
510 
511  float lxmatch = 9999.0;
512  float lymatch = 9999.0;
513  if (!expTrajMeasurements.empty()) {
514  if (glmatch < 9999.) { // if there is any propagated trajectory for this hit
515  const DetId &matchhit_detId = expTrajMeasurements[imatch].recHit()->geographicalId();
516 
517  int matchhit_ladder = PXBDetId(matchhit_detId).ladder();
518  int dladder = abs(matchhit_ladder - hit_ladder);
519  if (dladder > 10)
520  dladder = 20 - dladder;
521  LocalPoint lp = expTrajMeasurements[imatch].updatedState().localPosition();
522  lxmatch = fabs(lp.x() - chklp.x());
523  lymatch = fabs(lp.y() - chklp.y());
524  }
525  if (lxmatch < maxlxmatch_ && lymatch < maxlymatch_) {
526  if (testhit->getType() != TrackingRecHit::missing || keepOriginalMissingHit_) {
527  expTrajMeasurements.erase(expTrajMeasurements.begin() + imatch);
528  }
529  }
530 
531  } // expected trajectory measurment not empty
532  }
533  } // loop on trajectory measurments tmeasColl
534 
535  // if an extrapolated hit was found but not matched to an exisitng L1 hit
536  // then push the hit back into the collection now keep the first one that is
537  // left
538  if (!expTrajMeasurements.empty()) {
539  for (size_t f = 0; f < expTrajMeasurements.size(); f++) {
540  TrajectoryMeasurement AddHit = expTrajMeasurements[f];
541  if (AddHit.recHit()->getType() == TrackingRecHit::missing) {
542  tmeasColl.push_back(AddHit);
543  isBpixtrack = true;
544  }
545  }
546  }
547 
548  if (isBpixtrack || isFpixtrack) {
549  if (trackref->pt() < 0.6 || nStripHits < 8 || fabs(trackref->dxy(bestVtx->position())) > 0.05 ||
550  fabs(trackref->dz(bestVtx->position())) > 0.5)
551  continue;
552 
553  if (debug_) {
554  std::cout << "isBpixtrack : " << isBpixtrack << std::endl;
555  std::cout << "isFpixtrack : " << isFpixtrack << std::endl;
556  }
557  // std::cout<<"This tracks has so many hits:
558  // "<<tmeasColl.size()<<std::endl;
559  for (std::vector<TrajectoryMeasurement>::const_iterator tmeasIt = tmeasColl.begin(); tmeasIt != tmeasColl.end();
560  tmeasIt++) {
561  // if(! tmeasIt->updatedState().isValid()) continue;
562  TrajectoryStateOnSurface tsos = tmeasIt->updatedState();
563 
565  if (hit->geographicalId().det() != DetId::Tracker)
566  continue;
567  else {
568  // //residual
569  const DetId &hit_detId = hit->geographicalId();
570  // uint IntRawDetID = (hit_detId.rawId());
571  uint IntSubDetID = (hit_detId.subdetId());
572 
573  if (IntSubDetID == 0) {
574  if (debug_)
575  std::cout << "NO IntSubDetID\n";
576  continue;
577  }
578  if (IntSubDetID != PixelSubdetector::PixelBarrel && IntSubDetID != PixelSubdetector::PixelEndcap)
579  continue;
580 
581  int disk = 0;
582  int layer = 0;
583  int panel = 0;
584  int module = 0;
585  bool isHalfModule = false;
586 
587  const SiPixelRecHit *pixhit = dynamic_cast<const SiPixelRecHit *>(hit->hit());
588 
589  if (IntSubDetID == PixelSubdetector::PixelBarrel) { // it's a BPIX hit
590  layer = PixelBarrelName(hit_detId, pTT, isUpgrade).layerName();
591  isHalfModule = PixelBarrelName(hit_detId, pTT, isUpgrade).isHalfModule();
592 
593  if (hit->isValid()) { // fill the cluster probability in barrel
594  bool plus = true;
595  if ((PixelBarrelName(hit_detId, pTT, isUpgrade).shell() == PixelBarrelName::Shell::mO) ||
596  (PixelBarrelName(hit_detId, pTT, isUpgrade).shell() == PixelBarrelName::Shell::mI))
597  plus = false;
598  double clusterProbability = pixhit->clusterProbability(0);
599  if (clusterProbability != 0)
600  fillClusterProbability(layer, 0, plus, log10(clusterProbability));
601  }
602 
603  } else if (IntSubDetID == PixelSubdetector::PixelEndcap) { // it's an FPIX hit
604  disk = PixelEndcapName(hit_detId, pTT, isUpgrade).diskName();
605  panel = PixelEndcapName(hit_detId, pTT, isUpgrade).pannelName();
606  module = PixelEndcapName(hit_detId, pTT, isUpgrade).plaquetteName();
607 
608  if (hit->isValid()) {
609  bool plus = true;
610  if ((PixelEndcapName(hit_detId, pTT, isUpgrade).halfCylinder() == PixelEndcapName::HalfCylinder::mO) ||
611  (PixelEndcapName(hit_detId, pTT, isUpgrade).halfCylinder() == PixelEndcapName::HalfCylinder::mI))
612  plus = false;
613  double clusterProbability = pixhit->clusterProbability(0);
614  if (clusterProbability != 0)
615  fillClusterProbability(0, disk, plus, log10(clusterProbability));
616  }
617  }
618 
619  if (layer == 1) {
620  if (fabs(trackref->dxy(bestVtx->position())) > 0.01 || fabs(trackref->dz(bestVtx->position())) > 0.1)
621  continue;
622  if (!(L2hits > 0 && L3hits > 0) && !(L2hits > 0 && D1hits > 0) && !(D1hits > 0 && D2hits > 0))
623  continue;
624  } else if (layer == 2) {
625  if (fabs(trackref->dxy(bestVtx->position())) > 0.02 || fabs(trackref->dz(bestVtx->position())) > 0.1)
626  continue;
627  if (!(L1hits > 0 && L3hits > 0) && !(L1hits > 0 && D1hits > 0))
628  continue;
629  } else if (layer == 3) {
630  if (fabs(trackref->dxy(bestVtx->position())) > 0.02 || fabs(trackref->dz(bestVtx->position())) > 0.1)
631  continue;
632  if (!(L1hits > 0 && L2hits > 0))
633  continue;
634  } else if (layer == 4) {
635  if (fabs(trackref->dxy(bestVtx->position())) > 0.02 || fabs(trackref->dz(bestVtx->position())) > 0.1)
636  continue;
637  } else if (disk == 1) {
638  if (fabs(trackref->dxy(bestVtx->position())) > 0.05 || fabs(trackref->dz(bestVtx->position())) > 0.5)
639  continue;
640  if (!(L1hits > 0 && D2hits > 0) && !(L2hits > 0 && D2hits > 0))
641  continue;
642  } else if (disk == 2) {
643  if (fabs(trackref->dxy(bestVtx->position())) > 0.05 || fabs(trackref->dz(bestVtx->position())) > 0.5)
644  continue;
645  if (!(L1hits > 0 && D1hits > 0))
646  continue;
647  } else if (disk == 3) {
648  if (fabs(trackref->dxy(bestVtx->position())) > 0.05 || fabs(trackref->dz(bestVtx->position())) > 0.5)
649  continue;
650  }
651 
652  // check whether hit is valid or missing using track algo flag
653  bool isHitValid = hit->hit()->getType() == TrackingRecHit::valid;
654  bool isHitMissing = hit->hit()->getType() == TrackingRecHit::missing;
655 
656  if (debug_)
657  std::cout << "the hit is persistent\n";
658 
659  std::map<uint32_t, SiPixelHitEfficiencyModule *>::iterator pxd =
660  theSiPixelStructure.find((*hit).geographicalId().rawId());
661 
662  // calculate alpha and beta from cluster position
664  // LocalVector localDir = ltp.momentum()/ltp.momentum().mag();
665 
666  //*************** Edge cut ********************
667  double lx = tsos.localPosition().x();
668  double ly = tsos.localPosition().y();
669 
670  if (fabs(lx) > 0.55 || fabs(ly) > 3.0)
671  continue;
672 
673  bool passedFiducial = true;
674 
675  // Module fiducials:
676  if (IntSubDetID == PixelSubdetector::PixelBarrel && fabs(ly) >= 3.1)
677  passedFiducial = false;
678  if (IntSubDetID == PixelSubdetector::PixelEndcap &&
679  !((panel == 1 &&
680  ((module == 1 && fabs(ly) < 0.7) ||
681  ((module == 2 && fabs(ly) < 1.1) && !(disk == -1 && ly > 0.8 && lx > 0.2) &&
682  !(disk == 1 && ly < -0.7 && lx > 0.2) && !(disk == 2 && ly < -0.8)) ||
683  ((module == 3 && fabs(ly) < 1.5) && !(disk == -2 && lx > 0.1 && ly > 1.0) &&
684  !(disk == 2 && lx > 0.1 && ly < -1.0)) ||
685  ((module == 4 && fabs(ly) < 1.9) && !(disk == -2 && ly > 1.5) && !(disk == 2 && ly < -1.5)))) ||
686  (panel == 2 &&
687  ((module == 1 && fabs(ly) < 0.7) ||
688  (module == 2 && fabs(ly) < 1.2 && !(disk > 0 && ly > 1.1) && !(disk < 0 && ly < -1.1)) ||
689  (module == 3 && fabs(ly) < 1.6 && !(disk > 0 && ly > 1.5) && !(disk < 0 && ly < -1.5))))))
690  passedFiducial = false;
691  if (IntSubDetID == PixelSubdetector::PixelEndcap &&
692  ((panel == 1 && (module == 1 || (module >= 3 && abs(disk) == 1))) ||
693  (panel == 2 && ((module == 1 && abs(disk) == 2) || (module == 3 && abs(disk) == 1)))))
694  passedFiducial = false;
695  // ROC fiducials:
696  double ly_mod = fabs(ly);
697  if (IntSubDetID == PixelSubdetector::PixelEndcap && (panel + module) % 2 == 1)
698  ly_mod = ly_mod + 0.405;
699  float d_rocedge = fabs(fmod(ly_mod, 0.81) - 0.405);
700  if (d_rocedge <= 0.0625)
701  passedFiducial = false;
702  if (!((IntSubDetID == PixelSubdetector::PixelBarrel &&
703  ((!isHalfModule && fabs(lx) < 0.6) || (isHalfModule && lx > -0.3 && lx < 0.2))) ||
704  (IntSubDetID == PixelSubdetector::PixelEndcap &&
705  ((panel == 1 &&
706  ((module == 1 && fabs(lx) < 0.2) ||
707  (module == 2 && ((fabs(lx) < 0.55 && abs(disk) == 1) || (lx > -0.5 && lx < 0.2 && disk == -2) ||
708  (lx > -0.5 && lx < 0.0 && disk == 2))) ||
709  (module == 3 && lx > -0.6 && lx < 0.5) || (module == 4 && lx > -0.3 && lx < 0.15))) ||
710  (panel == 2 && ((module == 1 && fabs(lx) < 0.6) ||
711  (module == 2 && ((fabs(lx) < 0.55 && abs(disk) == 1) ||
712  (lx > -0.6 && lx < 0.5 && abs(disk) == 2))) ||
713  (module == 3 && fabs(lx) < 0.5)))))))
714  passedFiducial = false;
715  if (((IntSubDetID == PixelSubdetector::PixelBarrel && !isHalfModule) ||
716  (IntSubDetID == PixelSubdetector::PixelEndcap && !(panel == 1 && (module == 1 || module == 4)))) &&
717  fabs(lx) < 0.06)
718  passedFiducial = false;
719 
720  //*************** find closest clusters ********************
721  float dx_cl[2];
722  float dy_cl[2];
723  dx_cl[0] = dx_cl[1] = dy_cl[0] = dy_cl[1] = -9999.;
726  if (cpEstimator.isValid()) {
727  const PixelClusterParameterEstimator &cpe(*cpEstimator);
729  if (tracker.isValid()) {
730  const TrackerGeometry *tkgeom = &(*tracker);
731  edm::Handle<edmNew::DetSetVector<SiPixelCluster>> clusterCollectionHandle;
732  iEvent.getByToken(clusterCollectionToken_, clusterCollectionHandle);
733  if (clusterCollectionHandle.isValid()) {
734  const edmNew::DetSetVector<SiPixelCluster> &clusterCollection = *clusterCollectionHandle;
736  float minD[2];
737  minD[0] = minD[1] = 10000.;
738  for (; itClusterSet != clusterCollection.end(); itClusterSet++) {
739  DetId detId(itClusterSet->id());
740  if (detId.rawId() != hit->geographicalId().rawId())
741  continue;
742  // unsigned int sdId=detId.subdetId();
743  const PixelGeomDetUnit *pixdet = (const PixelGeomDetUnit *)tkgeom->idToDetUnit(detId);
744  edmNew::DetSet<SiPixelCluster>::const_iterator itCluster = itClusterSet->begin();
745  for (; itCluster != itClusterSet->end(); ++itCluster) {
746  LocalPoint lp(itCluster->x(), itCluster->y(), 0.);
748  lp = std::get<0>(params);
749  float D = sqrt((lp.x() - lx) * (lp.x() - lx) + (lp.y() - ly) * (lp.y() - ly));
750  if (D < minD[0]) {
751  minD[1] = minD[0];
752  dx_cl[1] = dx_cl[0];
753  dy_cl[1] = dy_cl[0];
754  minD[0] = D;
755  dx_cl[0] = lp.x();
756  dy_cl[0] = lp.y();
757  } else if (D < minD[1]) {
758  minD[1] = D;
759  dx_cl[1] = lp.x();
760  dy_cl[1] = lp.y();
761  }
762  } // loop on clusterSets
763  } // loop on clusterCollection
764  for (size_t i = 0; i < 2; i++) {
765  if (minD[i] < 9999.) {
766  dx_cl[i] = fabs(dx_cl[i] - lx);
767  dy_cl[i] = fabs(dy_cl[i] - ly);
768  }
769  }
770  } // valid clusterCollectionHandle
771  } // valid tracker
772  } // valid cpEstimator
773  // distance of hit from closest cluster!
774  float d_cl[2];
775  d_cl[0] = d_cl[1] = -9999.;
776  if (dx_cl[0] != -9999. && dy_cl[0] != -9999.)
777  d_cl[0] = sqrt(dx_cl[0] * dx_cl[0] + dy_cl[0] * dy_cl[0]);
778  if (dx_cl[1] != -9999. && dy_cl[1] != -9999.)
779  d_cl[1] = sqrt(dx_cl[1] * dx_cl[1] + dy_cl[1] * dy_cl[1]);
780  if (isHitMissing && (d_cl[0] < 0.05 || d_cl[1] < 0.05)) {
781  isHitMissing = false;
782  isHitValid = true;
783  }
784 
785  if (debug_) {
786  std::cout << "Ready to add hit in histogram:\n";
787  // std::cout << "detid: "<<hit_detId<<std::endl;
788  std::cout << "isHitValid: " << isHitValid << std::endl;
789  std::cout << "isHitMissing: " << isHitMissing << std::endl;
790  // std::cout << "passedEdgeCut: "<<passedFiducial<<std::endl;
791  }
792 
793  if (nStripHits < 11)
794  continue; // Efficiency plots are filled with hits on tracks that
795  // have at least 11 Strip hits
796 
797  if (pxd != theSiPixelStructure.end() && isHitValid && passedFiducial)
798  ++nvalid;
799  if (pxd != theSiPixelStructure.end() && isHitMissing && passedFiducial)
800  ++nmissing;
801 
802  if (pxd != theSiPixelStructure.end() && passedFiducial && (isHitValid || isHitMissing))
803  (*pxd).second->fill(pTT, ltp, isHitValid, modOn, ladOn, layOn, phiOn, bladeOn, diskOn, ringOn);
804 
805  //}//end if (persistent hit exists and is pixel hit)
806 
807  } // end of else
808 
809  } // end for (all traj measurements of pixeltrack)
810  } // end if (is pixeltrack)
811  else if (debug_)
812  std::cout << "no pixeltrack:\n";
813 
814  } // end loop on map entries
815 }
816 
817 DEFINE_FWK_MODULE(SiPixelHitEfficiencySource); // define this as a plug-in
edm::EDGetTokenT< MeasurementTrackerEvent > measurementTrackerEventToken_
Log< level::Info, true > LogVerbatim
uint32_t second(const DetId &id) const
static constexpr auto TEC
std::map< uint32_t, SiPixelHitEfficiencyModule * > theSiPixelStructure
virtual void setPropagationDirection(PropagationDirection dir)
Definition: Propagator.h:130
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
edm::ESGetToken< MeasurementTracker, CkfComponentsRecord > measurementTrackerToken_
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > trackerGeomTokenBeginRun_
SiPixelHitEfficiencySource(const edm::ParameterSet &)
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
const DetContainer & detsPXB() const
T z() const
Definition: PV3DBase.h:61
unsigned int module() const
det id
Definition: PXBDetId.h:37
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > clusterCollectionToken_
T const * product() const
Definition: Handle.h:70
virtual Propagator * clone() const =0
const DetContainer & detsPXF() const
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > trackerGeomToken_
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > trackerTopoTokenBeginRun_
const LocalTrajectoryParameters & localParameters() const
const DetContainer & dets() const override
Returm a vector of all GeomDet (including all GeomDetUnits)
data_type const * const_iterator
Definition: DetSetNew.h:31
virtual void fillClusterProbability(int, int, bool, double)
key_type key() const
Accessor for product key.
Definition: Ref.h:250
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
edm::EDGetTokenT< reco::VertexCollection > vertexCollectionToken_
const_iterator end(bool update=false) const
shell
Create symlink for executable/python cms config if needed.
virtual ReturnType getParameters(const SiPixelCluster &cl, const GeomDetUnit &det) const =0
const_iterator end() const
last iterator over the map (read only)
T getUntrackedParameter(std::string const &, T const &) const
void Fill(long long x)
size_type size() const
map size
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
int iEvent
Definition: GenABIO.cc:224
T const * product() const
Definition: ESHandle.h:86
GlobalPoint globalPosition() const
bool isHalfModule() const
full or half module
int diskName() const
disk id
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
T sqrt(T t)
Definition: SSEVec.h:19
bool setModuleFolder(const uint32_t &rawdetid=0, int type=0, bool isUpgrade=false)
Set folder name for a module or plaquette.
edm::EDGetTokenT< TrajTrackAssociationCollection > tracksrc_
HalfCylinder halfCylinder() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
id_type id(size_t cell) const
double f[11][100]
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
static constexpr auto TOB
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > trackerTopoToken_
float clusterProbability(unsigned int flags=0) const
Definition: SiPixelRecHit.cc:9
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
bool isValid() const
Definition: ESHandle.h:44
Shell shell() const
Log< level::Info, false > LogInfo
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
Definition: DetId.h:17
void analyze(const edm::Event &, const edm::EventSetup &) override
static constexpr auto TIB
const_iterator begin(bool update=false) const
DecomposeProduct< arg, typename Div::arg > D
Definition: Factorize.h:141
int layerName() const
layer id
bool isValid() const
Definition: HandleBase.h:70
edm::ESGetToken< Chi2MeasurementEstimatorBase, TrackingComponentsRecord > chi2MeasurementEstimatorBaseToken_
const_iterator begin() const
first iterator over the map (read only)
HLT enums.
edm::ESGetToken< Propagator, TrackingComponentsRecord > propagatorToken_
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
std::vector< TrajectoryMeasurement > measurements(const DetLayer &layer, const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est) const
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
modOn
online/offline RawDataErrors
unsigned int ladder() const
ladder id
Definition: PXBDetId.h:34
std::tuple< LocalPoint, LocalError, SiPixelRecHitQuality::QualWordType > ReturnType
static constexpr auto TID
void dqmBeginRun(const edm::Run &r, edm::EventSetup const &iSetup) override
int pannelName() const
pannel id
Definition: Run.h:45
ConstRecHitPointer const & recHit() const
Our base class.
Definition: SiPixelRecHit.h:23
virtual void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
edm::ESGetToken< PixelClusterParameterEstimator, TkPixelCPERecord > pixelClusterParameterEstimatorToken_
int plaquetteName() const
plaquetteId (in pannel)