CMS 3D CMS Logo

SiStripHitEfficiencyWorker.cc
Go to the documentation of this file.
1 // Package: CalibTracker/SiStripHitEfficiency
3 // Class: SiStripHitEfficiencyWorker
4 // Original Author: Pieter David
5 //
6 // Adapted from HitEff (Keith Ulmer -- University of Colorado, keith.ulmer@colorado.edu
7 // SiStripHitEffFromCalibTree (Christopher Edelmaier)
9 
63 
65 public:
66  explicit SiStripHitEfficiencyWorker(const edm::ParameterSet& conf);
67  ~SiStripHitEfficiencyWorker() override = default;
68  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
69 
70 private:
71  void bookHistograms(DQMStore::IBooker& booker, const edm::Run& run, const edm::EventSetup& setup) override;
72  void analyze(const edm::Event& e, const edm::EventSetup& c) override;
73  void fillForTraj(const TrajectoryAtInvalidHit& tm,
74  const TrackerTopology* tTopo,
75  const TrackerGeometry* tkgeom,
76  const StripClusterParameterEstimator& stripCPE,
77  const SiStripQuality& stripQuality,
78  const DetIdVector& fedErrorIds,
79  const edm::Handle<edm::DetSetVector<SiStripRawDigi>>& commonModeDigis,
80  const edmNew::DetSetVector<SiStripCluster>& theClusters,
81  int bunchCrossing,
82  float instLumi,
83  float PU,
84  bool highPurity);
85 
86  // ----------member data ---------------------------
88 
89  // event data tokens
100 
101  // event setup tokens
111 
112  // configurable parameters
114  unsigned int layers_;
115  bool DEBUG_;
116  bool addLumi_;
119  unsigned int trackMultiplicityCut_;
125  float resXSig_;
129  int bunchX_;
132  unsigned int nTEClayers_;
133 
134  // output file
135  std::set<uint32_t> badModules_;
136 
137  // for the missing hits recovery
138  std::vector<unsigned int> hitRecoveryCounters;
139  std::vector<unsigned int> hitTotalCounters;
141  std::vector<int> missHitPerLayer;
142 
143  struct EffME1 {
144  EffME1() : hTotal(nullptr), hFound(nullptr) {}
146 
147  void fill(double x, bool found, float weight = 1.) {
148  hTotal->Fill(x, weight);
149  if (found) {
150  hFound->Fill(x, weight);
151  }
152  }
153 
155  };
156  struct EffTkMap {
157  EffTkMap() : hTotal(nullptr), hFound(nullptr) {}
158  EffTkMap(std::unique_ptr<TkHistoMap>&& total, std::unique_ptr<TkHistoMap>&& found)
159  : hTotal(std::move(total)), hFound(std::move(found)) {}
160 
161  void fill(uint32_t id, bool found, float weight = 1.) {
162  hTotal->fill(id, weight);
163  if (found) {
164  hFound->fill(id, weight);
165  }
166  }
167 
168  bool check(uint32_t id) {
169  if (hTotal->getValue(id) < hFound->getValue(id)) {
170  return false;
171  } else {
172  return true;
173  }
174  }
175 
176  std::unique_ptr<TkHistoMap> hTotal, hFound;
177  };
178 
184  std::vector<MonitorElement*> h_resolution;
185  std::vector<EffME1> h_layer_vsLumi;
186  std::vector<EffME1> h_layer_vsBx;
187  std::vector<EffME1> h_layer_vsPU;
188  std::vector<EffME1> h_layer_vsCM;
189  std::vector<MonitorElement*> h_hotcold;
190 
192 };
193 
194 //
195 // constructors and destructor
196 //
197 
199  : scalerToken_(consumes<LumiScalersCollection>(conf.getParameter<edm::InputTag>("lumiScalers"))),
200  metaDataToken_(consumes<OnlineLuminosityRecord>(conf.getParameter<edm::InputTag>("metadata"))),
201  commonModeToken_(mayConsume<edm::DetSetVector<SiStripRawDigi>>(conf.getParameter<edm::InputTag>("commonMode"))),
202  combinatorialTracks_token_(
203  consumes<reco::TrackCollection>(conf.getParameter<edm::InputTag>("combinatorialTracks"))),
204  trajectories_token_(consumes<std::vector<Trajectory>>(conf.getParameter<edm::InputTag>("trajectories"))),
205  trajTrackAsso_token_(consumes<TrajTrackAssociationCollection>(conf.getParameter<edm::InputTag>("trajectories"))),
206  clusters_token_(
207  consumes<edmNew::DetSetVector<SiStripCluster>>(conf.getParameter<edm::InputTag>("siStripClusters"))),
208  digisCol_token_(consumes(conf.getParameter<edm::InputTag>("siStripDigis"))),
209  digisVec_token_(consumes(conf.getParameter<edm::InputTag>("siStripDigis"))),
210  trackerEvent_token_(consumes<MeasurementTrackerEvent>(conf.getParameter<edm::InputTag>("trackerEvent"))),
211  tTopoToken_(esConsumes()),
212  tkGeomToken_(esConsumes()),
213  stripCPEToken_(esConsumes(edm::ESInputTag{"", "StripCPEfromTrackAngle"})),
214  stripQualityToken_(esConsumes()),
215  magFieldToken_(esConsumes()),
216  measTrackerToken_(esConsumes()),
217  chi2EstimatorToken_(esConsumes(edm::ESInputTag{"", "Chi2"})),
218  propagatorToken_(esConsumes(edm::ESInputTag{"", "PropagatorWithMaterial"})),
219  tkDetMapToken_(esConsumes<edm::Transition::BeginRun>()),
220  dqmDir_(conf.getParameter<std::string>("dqmDir")),
221  layers_(conf.getParameter<int>("Layer")),
222  DEBUG_(conf.getUntrackedParameter<bool>("Debug", false)),
223  addLumi_(conf.getUntrackedParameter<bool>("addLumi", false)),
224  addCommonMode_(conf.getUntrackedParameter<bool>("addCommonMode", false)),
225  cutOnTracks_(conf.getParameter<bool>("cutOnTracks")),
226  trackMultiplicityCut_(conf.getParameter<unsigned int>("trackMultiplicity")),
227  useFirstMeas_(conf.getParameter<bool>("useFirstMeas")),
228  useLastMeas_(conf.getParameter<bool>("useLastMeas")),
229  useAllHitsFromTracksWithMissingHits_(conf.getParameter<bool>("useAllHitsFromTracksWithMissingHits")),
230  doMissingHitsRecovery_(conf.getParameter<bool>("doMissingHitsRecovery")),
231  clusterMatchingMethod_(conf.getParameter<int>("ClusterMatchingMethod")),
232  resXSig_(conf.getParameter<double>("ResXSig")),
233  clusterTracjDist_(conf.getParameter<double>("ClusterTrajDist")),
234  stripsApvEdge_(conf.getParameter<double>("StripsApvEdge")),
235  useOnlyHighPurityTracks_(conf.getParameter<bool>("UseOnlyHighPurityTracks")),
236  bunchX_(conf.getUntrackedParameter<int>("BunchCrossing", 0)),
237  showRings_(conf.getUntrackedParameter<bool>("ShowRings", false)),
238  showTOB6TEC9_(conf.getUntrackedParameter<bool>("ShowTOB6TEC9", false)) {
239  hitRecoveryCounters.resize(k_END_OF_LAYERS, 0);
240  hitTotalCounters.resize(k_END_OF_LAYERS, 0);
241  missHitPerLayer.resize(k_END_OF_LAYERS, 0);
242  totalNbHits = 0;
243 
244  nTEClayers_ = (showRings_ ? 7 : 9); // number of rings or wheels
245 
246  const std::string badModulesFile = conf.getUntrackedParameter<std::string>("BadModulesFile", "");
247  if (!badModulesFile.empty()) {
248  std::ifstream badModules_file(badModulesFile);
249  uint32_t badmodule_detid;
250  int mods, fiber1, fiber2, fiber3;
251  if (badModules_file.is_open()) {
253  while (getline(badModules_file, line)) {
254  if (badModules_file.eof())
255  continue;
256  std::stringstream ss(line);
257  ss >> badmodule_detid >> mods >> fiber1 >> fiber2 >> fiber3;
258  if (badmodule_detid != 0 && mods == 1 && (fiber1 == 1 || fiber2 == 1 || fiber3 == 1))
259  badModules_.insert(badmodule_detid);
260  }
261  badModules_file.close();
262  }
263  }
264  if (!badModules_.empty())
265  LogDebug("SiStripHitEfficiencyWorker") << "Remove additionnal bad modules from the analysis: ";
266  for (const auto badMod : badModules_) {
267  LogDebug("SiStripHitEfficiencyWorker") << " " << badMod;
268  }
269 }
270 
272  const edm::Run& run,
273  const edm::EventSetup& setup) {
274  booker.setCurrentFolder(fmt::format("{}/EventInfo", dqmDir_));
275  h_bx = booker.book1D("bx", "bx", 3600, 0, 3600);
276  h_instLumi = booker.book1D("instLumi", "inst. lumi.", 250, 0, 25000);
277  h_PU = booker.book1D("PU", "PU", 200, 0, 200);
278  h_nTracks = booker.book1D("ntracks", "n.tracks;n. tracks;n.events", 500, -0.5, 499.5);
279  h_nTracksVsPU = booker.bookProfile("nTracksVsPU", "n. tracks vs PU; PU; n.tracks ", 200, 0, 200, 500, -0.5, 499.5);
280 
281  calibData_.EventStats = booker.book2I("EventStats", "Statistics", 3, -0.5, 2.5, 1, 0, 1);
282  calibData_.EventStats->setBinLabel(1, "events count", 1);
283  calibData_.EventStats->setBinLabel(2, "tracks count", 1);
284  calibData_.EventStats->setBinLabel(3, "measurements count", 1);
285 
286  booker.setCurrentFolder(dqmDir_);
287  h_goodLayer = EffME1(booker.book1D("goodlayer_total", "goodlayer_total", 35, 0., 35.),
288  booker.book1D("goodlayer_found", "goodlayer_found", 35, 0., 35.));
289  h_allLayer = EffME1(booker.book1D("alllayer_total", "alllayer_total", 35, 0., 35.),
290  booker.book1D("alllayer_found", "alllayer_found", 35, 0., 35.));
291 
292  h_layer = EffME1(
293  booker.book1D(
294  "layer_found", "layer_found", bounds::k_END_OF_LAYERS, 0., static_cast<float>(bounds::k_END_OF_LAYERS)),
295  booker.book1D(
296  "layer_total", "layer_total", bounds::k_END_OF_LAYERS, 0., static_cast<float>(bounds::k_END_OF_LAYERS)));
297 
298  for (int layer = 1; layer != bounds::k_END_OF_LAYERS; ++layer) {
299  const auto lyrName = ::layerName(layer, showRings_, nTEClayers_);
300 
301  // book resolutions
302  booker.setCurrentFolder(fmt::format("{}/Resolutions", dqmDir_));
303  auto ihres = booker.book1D(Form("resol_layer_%i", layer), lyrName, 125, -125., 125.);
304  ihres->setAxisTitle("trajX-clusX [strip unit]");
305  h_resolution.push_back(ihres);
306 
307  // book plots vs Lumi
308  booker.setCurrentFolder(fmt::format("{}/VsLumi", dqmDir_));
309  h_layer_vsLumi.push_back(EffME1(booker.book1D(Form("layertotal_vsLumi_layer_%i", layer), lyrName, 100, 0, 25000),
310  booker.book1D(Form("layerfound_vsLumi_layer_%i", layer), lyrName, 100, 0, 25000)));
311 
312  // book plots vs Lumi
313  booker.setCurrentFolder(fmt::format("{}/VsPu", dqmDir_));
314  h_layer_vsPU.push_back(EffME1(booker.book1D(Form("layertotal_vsPU_layer_%i", layer), lyrName, 45, 0, 90),
315  booker.book1D(Form("layerfound_vsPU_layer_%i", layer), lyrName, 45, 0, 90)));
316  if (addCommonMode_) {
317  // book plots for common mode
318  booker.setCurrentFolder(fmt::format("{}/CommonMode", dqmDir_));
319  h_layer_vsCM.push_back(EffME1(booker.book1D(Form("layertotal_vsCM_layer_%i", layer), lyrName, 20, 0, 400),
320  booker.book1D(Form("layerfound_vsCM_layer_%i", layer), lyrName, 20, 0, 400)));
321  }
322 
323  // book plots vs Lumi
324  booker.setCurrentFolder(fmt::format("{}/VsBx", dqmDir_));
325  h_layer_vsBx.push_back(EffME1(
326  booker.book1D(Form("totalVsBx_layer%i", layer), Form("layer %i (%s)", layer, lyrName.c_str()), 3565, 0, 3565),
327  booker.book1D(Form("foundVsBx_layer%i", layer), Form("layer %i (%s)", layer, lyrName.c_str()), 3565, 0, 3565)));
328 
329  // book hot and cold
330  booker.setCurrentFolder(fmt::format("{}/MissingHits", dqmDir_));
331  if (layer <= bounds::k_LayersAtTOBEnd) {
332  const bool isTIB = layer <= bounds::k_LayersAtTIBEnd;
333  const auto partition = (isTIB ? "TIB" : "TOB");
334  const auto yMax = (isTIB ? 100 : 120);
335 
336  const auto& tit =
337  Form("%s%i: Map of missing hits", partition, (isTIB ? layer : layer - bounds::k_LayersAtTIBEnd));
338 
339  // histogram name must not contain ":" otherwise it fails upload to the GUI
340  // see https://github.com/cms-DQM/dqmgui_prod/blob/af0a388e8f57c60e51111585d298aeeea943367f/src/cpp/DQM/DQMStore.cc#L56
341  std::string name{tit};
342  ::replaceInString(name, ":", "");
343 
344  auto ihhotcold = booker.book2D(name, tit, 100, -1, 361, 100, -yMax, yMax);
345  ihhotcold->setAxisTitle("#phi [deg]", 1);
346  ihhotcold->setBinLabel(1, "360", 1);
347  ihhotcold->setBinLabel(50, "180", 1);
348  ihhotcold->setBinLabel(100, "0", 1);
349  ihhotcold->setAxisTitle("Global Z [cm]", 2);
350  ihhotcold->setOption("colz");
351  h_hotcold.push_back(ihhotcold);
352  } else {
353  const bool isTID = layer <= bounds::k_LayersAtTIDEnd;
354  const auto partitions =
355  (isTID ? std::vector<std::string>{"TIDplus", "TIDminus"} : std::vector<std::string>{"TECplus", "TECminus"});
356  const auto axMax = (isTID ? 100 : 120);
357  for (const auto& part : partitions) {
358  // create the title by replacing the minus/plus symbols
359  std::string forTitle{part};
360  ::replaceInString(forTitle, "minus", "-");
361  ::replaceInString(forTitle, "plus", "+");
362 
363  // histogram name must not contain ":" otherwise it fails upload to the GUI
364  // see https://github.com/cms-DQM/dqmgui_prod/blob/af0a388e8f57c60e51111585d298aeeea943367f/src/cpp/DQM/DQMStore.cc#L56
365  const auto& name = Form("%s %i Map of Missing Hits",
366  part.c_str(),
367  (isTID ? layer - bounds::k_LayersAtTOBEnd : layer - bounds::k_LayersAtTIDEnd));
368  const auto& tit = Form("%s%i: Map of Missing Hits",
369  forTitle.c_str(),
370  (isTID ? layer - bounds::k_LayersAtTOBEnd : layer - bounds::k_LayersAtTIDEnd));
371 
372  auto ihhotcold = booker.book2D(name, tit, 100, -axMax, axMax, 100, -axMax, axMax);
373  ihhotcold->setAxisTitle("Global Y", 1);
374  ihhotcold->setBinLabel(1, "+Y", 1);
375  ihhotcold->setBinLabel(50, "0", 1);
376  ihhotcold->setBinLabel(100, "-Y", 1);
377  ihhotcold->setAxisTitle("Global X", 2);
378  ihhotcold->setBinLabel(1, "-X", 2);
379  ihhotcold->setBinLabel(50, "0", 2);
380  ihhotcold->setBinLabel(100, "+X", 2);
381  ihhotcold->setOption("colz");
382  h_hotcold.push_back(ihhotcold);
383  }
384  }
385  }
386 
387  // come back to the main folder
388  booker.setCurrentFolder(dqmDir_);
389  const auto tkDetMapFolder = fmt::format("{}/TkDetMaps", dqmDir_);
390 
391  const TkDetMap* tkDetMap = &setup.getData(tkDetMapToken_);
392  h_module =
393  EffTkMap(std::make_unique<TkHistoMap>(tkDetMap, booker, tkDetMapFolder, "perModule_total", 0, false, true),
394  std::make_unique<TkHistoMap>(tkDetMap, booker, tkDetMapFolder, "perModule_found", 0, false, true));
395 
396  // fill the FED Errors
397  booker.setCurrentFolder(dqmDir_);
398  const auto FEDErrorMapFolder = fmt::format("{}/FEDErrorTkDetMaps", dqmDir_);
400  std::make_unique<TkHistoMap>(tkDetMap, booker, FEDErrorMapFolder, "perModule_FEDErrors", 0, false, true);
401 }
402 
404  const auto tTopo = &es.getData(tTopoToken_);
405 
406  // bool DEBUG_ = false;
407 
408  LogDebug("SiStripHitEfficiencyWorker") << "beginning analyze from HitEff";
409 
410  // Step A: Get Inputs
411 
412  // Luminosity informations
415 
416  float instLumi = 0;
417  float PU = 0;
418  if (addLumi_) {
419  if (lumiScalers.isValid() && !lumiScalers->empty()) {
420  if (lumiScalers->begin() != lumiScalers->end()) {
421  instLumi = lumiScalers->begin()->instantLumi();
422  PU = lumiScalers->begin()->pileup();
423  }
424  } else if (metaData.isValid()) {
425  instLumi = metaData->instLumi();
426  PU = metaData->avgPileUp();
427  } else {
428  edm::LogWarning("SiStripHitEfficiencyWorker") << "could not find a source for the Luminosity and PU";
429  }
430  }
431 
432  h_bx->Fill(e.bunchCrossing());
434  h_PU->Fill(PU);
435 
437  if (addCommonMode_)
438  e.getByToken(commonModeToken_, commonModeDigis);
439 
441  e.getByToken(combinatorialTracks_token_, tracksCKF);
442 
443  edm::Handle<std::vector<Trajectory>> TrajectoryCollectionCKF;
444  e.getByToken(trajectories_token_, TrajectoryCollectionCKF);
445 
446  edm::Handle<TrajTrackAssociationCollection> trajTrackAssociationHandle;
447  e.getByToken(trajTrackAsso_token_, trajTrackAssociationHandle);
448 
450  e.getByToken(clusters_token_, theClusters);
451 
452  // get the list of module IDs with FED-detected errors
453  // - In Aug-2023, the data format was changed from DetIdCollection to DetIdVector.
454  // - To provide some level of backward-compatibility,
455  // the plugin checks for both types giving preference to the new format.
456  // - If only the old format is available, the collection is
457  // converted to the new format, then used downstream.
458  auto const& fedErrorIdsCol_h = e.getHandle(digisCol_token_);
459  auto const& fedErrorIdsVec_h = e.getHandle(digisVec_token_);
460  if (not fedErrorIdsCol_h.isValid() and not fedErrorIdsVec_h.isValid()) {
461  throw cms::Exception("InvalidProductSiStripDetIdsWithFEDErrors")
462  << "no valid product for SiStrip DetIds with FED errors (see parameter \"siStripDigis\"), "
463  "neither for new format (DetIdVector) nor old format (DetIdCollection)";
464  }
465  auto const& fedErrorIds = fedErrorIdsVec_h.isValid() ? *fedErrorIdsVec_h : fedErrorIdsCol_h->as_vector();
466 
467  // fill the calibData with the FEDErrors
468  for (const auto& fedErr : fedErrorIds) {
469  // fill the TkHistoMap occupancy map
470  calibData_.FEDErrorOccupancy->fill(fedErr.rawId(), 1.);
471 
472  // fill the unordered map
473  if (calibData_.fedErrorCounts.find(fedErr.rawId()) != calibData_.fedErrorCounts.end()) {
474  calibData_.fedErrorCounts[fedErr.rawId()] += 1;
475  } else {
476  calibData_.fedErrorCounts.insert(std::make_pair(fedErr.rawId(), 1));
477  }
478  }
479 
482 
483  const auto tkgeom = &es.getData(tkGeomToken_);
484  const auto& stripcpe = es.getData(stripCPEToken_);
485  const auto& stripQuality = es.getData(stripQualityToken_);
486  const auto& magField = es.getData(magFieldToken_);
487  const auto& measTracker = es.getData(measTrackerToken_);
488  const auto& chi2Estimator = es.getData(chi2EstimatorToken_);
489  const auto& prop = es.getData(propagatorToken_);
490 
491  // Tracking
492  LogDebug("SiStripHitEfficiencyWorker") << "number ckf tracks found = " << tracksCKF->size();
493 
494  h_nTracks->Fill(tracksCKF->size());
495  h_nTracksVsPU->Fill(PU, tracksCKF->size());
496 
497  // bin 0: one entry for each event
498  calibData_.EventStats->Fill(0., 0., 1);
499  // bin 1: one entry for each track
500  calibData_.EventStats->Fill(1., 0., tracksCKF->size());
501 
502  if (!tracksCKF->empty()) {
503  if (cutOnTracks_ && (tracksCKF->size() >= trackMultiplicityCut_))
504  return;
505  if (cutOnTracks_)
506  LogDebug("SiStripHitEfficiencyWorker")
507  << "starting checking good event with < " << trackMultiplicityCut_ << " tracks";
508 
509  // actually should do a loop over all the tracks in the event here
510 
511  // Looping over traj-track associations to be able to get traj & track informations
512  for (const auto& trajTrack : *trajTrackAssociationHandle) {
513  // for each track, fill some variables such as number of hits and momentum
514 
515  const bool highPurity = trajTrack.val->quality(reco::TrackBase::TrackQuality::highPurity);
516  auto TMeas = trajTrack.key->measurements();
517  totalNbHits += int(TMeas.size());
518 
519  /*
520  const bool hasMissingHits = std::any_of(std::begin(TMeas), std::end(TMeas), [](const auto& tm) {
521  return tm.recHit()->getType() == TrackingRecHit::Type::missing;
522  });
523  */
524 
525  // Check whether the trajectory has some missing hits
526  bool hasMissingHits{false};
527  int previous_layer{999};
528  std::vector<unsigned int> missedLayers;
529 
530  for (const auto& itm : TMeas) {
531  auto theHit = itm.recHit();
532  unsigned int iidd = theHit->geographicalId().rawId();
533  int layer = ::checkLayer(iidd, tTopo);
534  int missedLayer = layer + 1;
535  int diffPreviousLayer = (layer - previous_layer);
537  //Layers from TIB + TOB
538  if (diffPreviousLayer == -2 && missedLayer > k_LayersStart && missedLayer < k_LayersAtTOBEnd) {
539  missHitPerLayer[missedLayer] += 1;
540  hasMissingHits = true;
541  }
542  //Layers from TID
543  else if (diffPreviousLayer == -2 && (missedLayer > k_LayersAtTOBEnd + 1 && missedLayer <= k_LayersAtTIDEnd)) {
544  missHitPerLayer[missedLayer] += 1;
545  hasMissingHits = true;
546  }
547  //Layers from TEC
548  else if (diffPreviousLayer == -2 && missedLayer > k_LayersAtTIDEnd && missedLayer <= k_LayersAtTECEnd) {
549  missHitPerLayer[missedLayer] += 1;
550  hasMissingHits = true;
551  }
552 
553  //##### TID Layer 11 in transition TID -> TIB : layer is in TIB, previous layer = 12
554  if ((layer > k_LayersStart && layer <= k_LayersAtTIBEnd) && (previous_layer == 12)) {
555  missHitPerLayer[11] += 1;
556  hasMissingHits = true;
557  }
558 
559  //##### TEC Layer 14 in transition TEC -> TOB : layer is in TOB, previous layer = 15
560  if ((layer > k_LayersAtTIBEnd && layer <= k_LayersAtTOBEnd) && (previous_layer == 15)) {
561  missHitPerLayer[14] += 1;
562  hasMissingHits = true;
563  }
564  }
565  if (theHit->getType() == TrackingRecHit::Type::missing)
566  hasMissingHits = true;
567 
568  if (hasMissingHits)
569  missedLayers.push_back(layer);
570  previous_layer = layer;
571  }
572 
573  // Loop on each measurement and take into consideration
574  //--------------------------------------------------------
575  unsigned int prev_TKlayers = 0;
576  for (auto itm = TMeas.cbegin(); itm != TMeas.cend(); ++itm) {
577  const auto theInHit = (*itm).recHit();
578 
579  //bin 2: one entry for each measurement
580  calibData_.EventStats->Fill(2., 0., 1.);
581 
582  LogDebug("SiStripHitEfficiencyWorker") << "theInHit is valid = " << theInHit->isValid();
583 
584  unsigned int iidd = theInHit->geographicalId().rawId();
585 
586  unsigned int TKlayers = ::checkLayer(iidd, tTopo);
587 
588  // do not bother with pixel hits
589  if (DetId(iidd).subdetId() < SiStripSubdetector::TIB)
590  continue;
591 
592  LogDebug("SiStripHitEfficiencyWorker") << "TKlayer from trajectory: " << TKlayers << " from module = " << iidd
593  << " matched/stereo/rphi = " << ((iidd & 0x3) == 0) << "/"
594  << ((iidd & 0x3) == 1) << "/" << ((iidd & 0x3) == 2);
595 
596  // Test first and last points of the trajectory
597  // the list of measurements starts from outer layers !!! This could change -> should add a check
598  bool isFirstMeas = (itm == (TMeas.end() - 1));
599  bool isLastMeas = (itm == (TMeas.begin()));
600 
601  if (!useFirstMeas_ && isFirstMeas)
602  continue;
603  if (!useLastMeas_ && isLastMeas)
604  continue;
605 
606  // In case of missing hit in the track, check whether to use the other hits or not.
607  if (hasMissingHits && theInHit->getType() != TrackingRecHit::Type::missing &&
609  continue;
610 
611  // If Trajectory measurement from TOB 6 or TEC 9, skip it because it's always valid they are filled later
612  if (TKlayers == bounds::k_LayersAtTOBEnd || TKlayers == bounds::k_LayersAtTECEnd) {
613  LogDebug("SiStripHitEfficiencyWorker") << "skipping original TM for TOB 6 or TEC 9";
614  continue;
615  }
616 
617  std::vector<TrajectoryAtInvalidHit> TMs;
618 
619  // Make AnalyticalPropagat // TODO where to save these?or to use in TAVH constructor
621 
622  // for double sided layers check both sensors--if no hit was found on either sensor surface,
623  // the trajectory measurements only have one invalid hit entry on the matched surface
624  // so get the TrajectoryAtInvalidHit for both surfaces and include them in the study
625  if (::isDoubleSided(iidd, tTopo) && ((iidd & 0x3) == 0)) {
626  // do hit eff check twice--once for each sensor
627  //add a TM for each surface
628  TMs.emplace_back(*itm, tTopo, tkgeom, propagator, 1);
629  TMs.emplace_back(*itm, tTopo, tkgeom, propagator, 2);
630  } else if (::isDoubleSided(iidd, tTopo) && (!::check2DPartner(iidd, TMeas))) {
631  // if only one hit was found the trajectory measurement is on that sensor surface, and the other surface from
632  // the matched layer should be added to the study as well
633  TMs.emplace_back(*itm, tTopo, tkgeom, propagator, 1);
634  TMs.emplace_back(*itm, tTopo, tkgeom, propagator, 2);
635  LogDebug("SiStripHitEfficiencyWorker") << " found a hit with a missing partner";
636  } else {
637  //only add one TM for the single surface and the other will be added in the next iteration
638  TMs.emplace_back(*itm, tTopo, tkgeom, propagator);
639  }
640 
641  bool missingHitAdded{false};
642  std::vector<TrajectoryMeasurement> tmpTmeas;
643  unsigned int misLayer = TKlayers + 1;
644  //Use bool doMissingHitsRecovery to add possible missing hits based on actual/previous hit
646  if (int(TKlayers) - int(prev_TKlayers) == -2) {
647  const DetLayer* detlayer = itm->layer();
648  const LayerMeasurements layerMeasurements{measTracker, *measurementTrackerEvent};
649  const TrajectoryStateOnSurface tsos = itm->updatedState();
650  std::vector<DetLayer::DetWithState> compatDets = detlayer->compatibleDets(tsos, prop, chi2Estimator);
651 
652  if (misLayer > k_LayersAtTIDEnd && misLayer < k_LayersAtTECEnd) { //TEC
653  std::vector<ForwardDetLayer const*> negTECLayers = measTracker.geometricSearchTracker()->negTecLayers();
654  std::vector<ForwardDetLayer const*> posTECLayers = measTracker.geometricSearchTracker()->posTecLayers();
655  const DetLayer* tecLayerneg = negTECLayers[misLayer - k_LayersAtTIDEnd - 1];
656  const DetLayer* tecLayerpos = posTECLayers[misLayer - k_LayersAtTIDEnd - 1];
657  if (tTopo->tecSide(iidd) == 1) {
658  tmpTmeas = layerMeasurements.measurements(*tecLayerneg, tsos, prop, chi2Estimator);
659  } else if (tTopo->tecSide(iidd) == 2) {
660  tmpTmeas = layerMeasurements.measurements(*tecLayerpos, tsos, prop, chi2Estimator);
661  }
662  }
663 
664  else if (misLayer == (k_LayersAtTIDEnd - 1) ||
665  misLayer == k_LayersAtTIDEnd) { // This is for TID layers 12 and 13
666  std::vector<ForwardDetLayer const*> negTIDLayers = measTracker.geometricSearchTracker()->negTidLayers();
667  std::vector<ForwardDetLayer const*> posTIDLayers = measTracker.geometricSearchTracker()->posTidLayers();
668  const DetLayer* tidLayerneg = negTIDLayers[misLayer - k_LayersAtTOBEnd - 1];
669  const DetLayer* tidLayerpos = posTIDLayers[misLayer - k_LayersAtTOBEnd - 1];
670 
671  if (tTopo->tidSide(iidd) == 1) {
672  tmpTmeas = layerMeasurements.measurements(*tidLayerneg, tsos, prop, chi2Estimator);
673  } else if (tTopo->tidSide(iidd) == 2) {
674  tmpTmeas = layerMeasurements.measurements(*tidLayerpos, tsos, prop, chi2Estimator);
675  }
676  }
677 
678  if (misLayer > k_LayersStart && misLayer < k_LayersAtTOBEnd) { // Barrel
679 
680  std::vector<BarrelDetLayer const*> barrelTIBLayers = measTracker.geometricSearchTracker()->tibLayers();
681  std::vector<BarrelDetLayer const*> barrelTOBLayers = measTracker.geometricSearchTracker()->tobLayers();
682 
683  if (misLayer > k_LayersStart && misLayer <= k_LayersAtTIBEnd) {
684  const DetLayer* tibLayer = barrelTIBLayers[misLayer - k_LayersStart - 1];
685  tmpTmeas = layerMeasurements.measurements(*tibLayer, tsos, prop, chi2Estimator);
686  } else if (misLayer > k_LayersAtTIBEnd && misLayer < k_LayersAtTOBEnd) {
687  const DetLayer* tobLayer = barrelTOBLayers[misLayer - k_LayersAtTIBEnd - 1];
688  tmpTmeas = layerMeasurements.measurements(*tobLayer, tsos, prop, chi2Estimator);
689  }
690  }
691  }
692  if ((int(TKlayers) > k_LayersStart && int(TKlayers) <= k_LayersAtTIBEnd) && int(prev_TKlayers) == 12) {
693  const DetLayer* detlayer = itm->layer();
694  const LayerMeasurements layerMeasurements{measTracker, *measurementTrackerEvent};
695  const TrajectoryStateOnSurface tsos = itm->updatedState();
696  std::vector<DetLayer::DetWithState> compatDets = detlayer->compatibleDets(tsos, prop, chi2Estimator);
697  std::vector<ForwardDetLayer const*> negTIDLayers = measTracker.geometricSearchTracker()->negTidLayers();
698  std::vector<ForwardDetLayer const*> posTIDLayers = measTracker.geometricSearchTracker()->posTidLayers();
699 
700  const DetLayer* tidLayerneg = negTIDLayers[k_LayersStart];
701  const DetLayer* tidLayerpos = posTIDLayers[k_LayersStart];
702  if (tTopo->tidSide(iidd) == 1) {
703  tmpTmeas = layerMeasurements.measurements(*tidLayerneg, tsos, prop, chi2Estimator);
704  } else if (tTopo->tidSide(iidd) == 2) {
705  tmpTmeas = layerMeasurements.measurements(*tidLayerpos, tsos, prop, chi2Estimator);
706  }
707  }
708 
709  if ((int(TKlayers) > k_LayersAtTIBEnd && int(TKlayers) <= k_LayersAtTOBEnd) && int(prev_TKlayers) == 15) {
710  const DetLayer* detlayer = itm->layer();
711  const LayerMeasurements layerMeasurements{measTracker, *measurementTrackerEvent};
712  const TrajectoryStateOnSurface tsos = itm->updatedState();
713  std::vector<DetLayer::DetWithState> compatDets = detlayer->compatibleDets(tsos, prop, chi2Estimator);
714 
715  std::vector<ForwardDetLayer const*> negTECLayers = measTracker.geometricSearchTracker()->negTecLayers();
716  std::vector<ForwardDetLayer const*> posTECLayers = measTracker.geometricSearchTracker()->posTecLayers();
717 
718  const DetLayer* tecLayerneg = negTECLayers[k_LayersStart];
719  const DetLayer* tecLayerpos = posTECLayers[k_LayersStart];
720  if (tTopo->tecSide(iidd) == 1) {
721  tmpTmeas = layerMeasurements.measurements(*tecLayerneg, tsos, prop, chi2Estimator);
722  } else if (tTopo->tecSide(iidd) == 2) {
723  tmpTmeas = layerMeasurements.measurements(*tecLayerpos, tsos, prop, chi2Estimator);
724  }
725  }
726 
727  if (!tmpTmeas.empty()) {
728  TrajectoryMeasurement TM_tmp(tmpTmeas.back());
729  unsigned int iidd_tmp = TM_tmp.recHit()->geographicalId().rawId();
730  if (iidd_tmp != 0) {
731  LogDebug("SiStripHitEfficiency:HitEff") << " hit actually being added to TM vector";
732  if ((!useAllHitsFromTracksWithMissingHits_ || (!useFirstMeas_ && isFirstMeas)))
733  TMs.clear();
734  if (::isDoubleSided(iidd_tmp, tTopo)) {
735  TMs.push_back(TrajectoryAtInvalidHit(TM_tmp, tTopo, tkgeom, propagator, 1));
736  TMs.push_back(TrajectoryAtInvalidHit(TM_tmp, tTopo, tkgeom, propagator, 2));
737  } else
738  TMs.push_back(TrajectoryAtInvalidHit(TM_tmp, tTopo, tkgeom, propagator));
739  missingHitAdded = true;
740  hitRecoveryCounters[misLayer] += 1;
741  }
742  }
743  }
744 
745  prev_TKlayers = TKlayers;
746  if (!useFirstMeas_ && isFirstMeas && !missingHitAdded)
747  continue;
748  if (!useLastMeas_ && isLastMeas)
749  continue;
750  bool hitsWithBias = false;
751  for (auto ilayer : missedLayers) {
752  if (ilayer < TKlayers)
753  hitsWithBias = true;
754  }
755  if (hasMissingHits && theInHit->getType() != TrackingRecHit::Type::missing && !missingHitAdded &&
756  hitsWithBias && !useAllHitsFromTracksWithMissingHits_) {
757  continue;
758  }
759 
761  //Now check for tracks at TOB6 and TEC9
762 
763  // to make sure we only propagate on the last TOB5 hit check the next entry isn't also in TOB5
764  // to avoid bias, make sure the TOB5 hit is valid (an invalid hit on TOB5 could only exist with a valid hit on TOB6)
765  const auto nextId = (itm + 1 != TMeas.end()) ? (itm + 1)->recHit()->geographicalId() : DetId{}; // null if last
766 
767  if (TKlayers == 9 && theInHit->isValid() && !((!nextId.null()) && (::checkLayer(nextId.rawId(), tTopo) == 9))) {
768  // if ( TKlayers==9 && itm==TMeas.rbegin()) {
769  // if ( TKlayers==9 && (itm==TMeas.back()) ) { // to check for only the last entry in the trajectory for propagation
770  const DetLayer* tob6 = measTracker.geometricSearchTracker()->tobLayers().back();
771  const LayerMeasurements theLayerMeasurements{measTracker, *measurementTrackerEvent};
772  const TrajectoryStateOnSurface tsosTOB5 = itm->updatedState();
773  const auto tmp = theLayerMeasurements.measurements(*tob6, tsosTOB5, prop, chi2Estimator);
774 
775  if (!tmp.empty()) {
776  LogDebug("SiStripHitEfficiencyWorker") << "size of TM from propagation = " << tmp.size();
777 
778  // take the last of the TMs, which is always an invalid hit
779  // if no detId is available, ie detId==0, then no compatible layer was crossed
780  // otherwise, use that TM for the efficiency measurement
781  const auto& tob6TM = tmp.back();
782  const auto& tob6Hit = tob6TM.recHit();
783  if (tob6Hit->geographicalId().rawId() != 0) {
784  LogDebug("SiStripHitEfficiencyWorker") << "tob6 hit actually being added to TM vector";
785  TMs.emplace_back(tob6TM, tTopo, tkgeom, propagator);
786  }
787  }
788  }
789 
790  // same for TEC8
791  if (TKlayers == 21 && theInHit->isValid() &&
792  !((!nextId.null()) && (::checkLayer(nextId.rawId(), tTopo) == 21))) {
793  const DetLayer* tec9pos = measTracker.geometricSearchTracker()->posTecLayers().back();
794  const DetLayer* tec9neg = measTracker.geometricSearchTracker()->negTecLayers().back();
795 
796  const LayerMeasurements theLayerMeasurements{measTracker, *measurementTrackerEvent};
797  const TrajectoryStateOnSurface tsosTEC9 = itm->updatedState();
798 
799  // check if track on positive or negative z
800  if (!(iidd == SiStripSubdetector::TEC))
801  LogDebug("SiStripHitEfficiencyWorker") << "there is a problem with TEC 9 extrapolation";
802 
803  //LogDebug("SiStripHitEfficiencyWorker") << " tec9 id = " << iidd << " and side = " << tTopo->tecSide(iidd) ;
804  std::vector<TrajectoryMeasurement> tmp;
805  if (tTopo->tecSide(iidd) == 1) {
806  tmp = theLayerMeasurements.measurements(*tec9neg, tsosTEC9, prop, chi2Estimator);
807  //LogDebug("SiStripHitEfficiencyWorker") << "on negative side" ;
808  }
809  if (tTopo->tecSide(iidd) == 2) {
810  tmp = theLayerMeasurements.measurements(*tec9pos, tsosTEC9, prop, chi2Estimator);
811  //LogDebug("SiStripHitEfficiencyWorker") << "on positive side" ;
812  }
813 
814  if (!tmp.empty()) {
815  // take the last of the TMs, which is always an invalid hit
816  // if no detId is available, ie detId==0, then no compatible layer was crossed
817  // otherwise, use that TM for the efficiency measurement
818  const auto& tec9TM = tmp.back();
819  const auto& tec9Hit = tec9TM.recHit();
820 
821  const unsigned int tec9id = tec9Hit->geographicalId().rawId();
822  LogDebug("SiStripHitEfficiencyWorker")
823  << "tec9id = " << tec9id << " is Double sided = " << ::isDoubleSided(tec9id, tTopo)
824  << " and 0x3 = " << (tec9id & 0x3);
825 
826  if (tec9Hit->geographicalId().rawId() != 0) {
827  LogDebug("SiStripHitEfficiencyWorker") << "tec9 hit actually being added to TM vector";
828  // in tec the hit can be single or doubled sided. whenever the invalid hit at the end of vector of TMs is
829  // double sided it is always on the matched surface, so we need to split it into the true sensor surfaces
830  if (::isDoubleSided(tec9id, tTopo)) {
831  TMs.emplace_back(tec9TM, tTopo, tkgeom, propagator, 1);
832  TMs.emplace_back(tec9TM, tTopo, tkgeom, propagator, 2);
833  } else
834  TMs.emplace_back(tec9TM, tTopo, tkgeom, propagator);
835  }
836  } //else LogDebug("SiStripHitEfficiencyWorker") << "tec9 tmp empty" ;
837  }
838 
839  for (const auto& tm : TMs) {
840  fillForTraj(tm,
841  tTopo,
842  tkgeom,
843  stripcpe,
844  stripQuality,
845  fedErrorIds,
846  commonModeDigis,
847  *theClusters,
848  e.bunchCrossing(),
849  instLumi,
850  PU,
851  highPurity);
852  }
853  LogDebug("SiStripHitEfficiencyWorker") << "After looping over TrajAtValidHit list";
854  }
855  LogDebug("SiStripHitEfficiencyWorker") << "end TMeasurement loop";
856  }
857  LogDebug("SiStripHitEfficiencyWorker") << "end of trajectories loop";
858  }
859 }
860 
862  const TrackerTopology* tTopo,
863  const TrackerGeometry* tkgeom,
864  const StripClusterParameterEstimator& stripCPE,
865  const SiStripQuality& stripQuality,
866  const DetIdVector& fedErrorIds,
867  const edm::Handle<edm::DetSetVector<SiStripRawDigi>>& commonModeDigis,
868  const edmNew::DetSetVector<SiStripCluster>& theClusters,
869  int bunchCrossing,
870  float instLumi,
871  float PU,
872  bool highPurity) {
873  // --> Get trajectory from combinatedStat& e
874  const auto iidd = tm.monodet_id();
875  LogDebug("SiStripHitEfficiencyWorker") << "setting iidd = " << iidd << " before checking efficiency and ";
876 
877  const auto xloc = tm.localX();
878  const auto yloc = tm.localY();
879 
880  const auto xErr = tm.localErrorX();
881  const auto yErr = tm.localErrorY();
882 
883  int TrajStrip = -1;
884 
885  // reget layer from iidd here, to account for TOB 6 and TEC 9 TKlayers being off
886  const auto TKlayers = ::checkLayer(iidd, tTopo);
887 
888  const bool withinAcceptance =
889  tm.withinAcceptance() && (!::isInBondingExclusionZone(iidd, TKlayers, yloc, yErr, tTopo));
890 
891  if ( // (TKlayers > 0) && // FIXME confirm this
892  ((layers_ == TKlayers) ||
893  (layers_ == bounds::k_LayersStart))) { // Look at the layer not used to reconstruct the track
894  LogDebug("SiStripHitEfficiencyWorker") << "Looking at layer under study";
895  unsigned int ModIsBad = 2;
896  unsigned int SiStripQualBad = 0;
897  float commonMode = -100;
898 
899  // RPhi RecHit Efficiency
900 
901  if (!theClusters.empty()) {
902  LogDebug("SiStripHitEfficiencyWorker") << "Checking clusters with size = " << theClusters.size();
903  std::vector<::ClusterInfo> VCluster_info; //fill with X residual, X residual pull, local X
904  const auto idsv = theClusters.find(iidd);
905  if (idsv != theClusters.end()) {
906  //if (DEBUG_) LogDebug("SiStripHitEfficiencyWorker") << "the ID from the dsv = " << dsv.id() ;
907  LogDebug("SiStripHitEfficiencyWorker")
908  << "found (ClusterId == iidd) with ClusterId = " << idsv->id() << " and iidd = " << iidd;
909  const auto stripdet = dynamic_cast<const StripGeomDetUnit*>(tkgeom->idToDetUnit(DetId(iidd)));
910  const StripTopology& Topo = stripdet->specificTopology();
911 
912  float hbedge = 0.0;
913  float htedge = 0.0;
914  float hapoth = 0.0;
915  float uylfac = 0.0;
916  float uxlden = 0.0;
917  if (TKlayers > bounds::k_LayersAtTOBEnd) {
918  const BoundPlane& plane = stripdet->surface();
919  const TrapezoidalPlaneBounds* trapezoidalBounds(
920  dynamic_cast<const TrapezoidalPlaneBounds*>(&(plane.bounds())));
921  std::array<const float, 4> const& parameterTrap = (*trapezoidalBounds).parameters(); // el bueno aqui
922  hbedge = parameterTrap[0];
923  htedge = parameterTrap[1];
924  hapoth = parameterTrap[3];
925  uylfac = (htedge - hbedge) / (htedge + hbedge) / hapoth;
926  uxlden = 1 + yloc * uylfac;
927  }
928 
929  // Need to know position of trajectory in strip number for selecting the right APV later
930  if (TrajStrip == -1) {
931  int nstrips = Topo.nstrips();
932  float pitch = stripdet->surface().bounds().width() / nstrips;
933  TrajStrip = xloc / pitch + nstrips / 2.0;
934  // Need additionnal corrections for endcap
935  if (TKlayers > bounds::k_LayersAtTOBEnd) {
936  const float TrajLocXMid = xloc / (1 + (htedge - hbedge) * yloc / (htedge + hbedge) /
937  hapoth); // radialy extrapolated x loc position at middle
938  TrajStrip = TrajLocXMid / pitch + nstrips / 2.0;
939  }
940  //LogDebug("SiStripHitEfficiency")<<" Layer "<<TKlayers<<" TrajStrip: "<<nstrips<<" "<<pitch<<" "<<TrajStrip;;
941  }
942 
943  for (const auto& clus : *idsv) {
945  float res = (parameters.first.x() - xloc);
946  float sigma = ::checkConsistency(parameters, xloc, xErr);
947  // The consistency is probably more accurately measured with the Chi2MeasurementEstimator. To use it
948  // you need a TransientTrackingRecHit instead of the cluster
949  //theEstimator= new Chi2MeasurementEstimator(30);
950  //const Chi2MeasurementEstimator *theEstimator(100);
951  //theEstimator->estimate(tm.tsos(), TransientTrackingRecHit);
952 
953  if (TKlayers > bounds::k_LayersAtTOBEnd) {
954  res = parameters.first.x() - xloc / uxlden; // radialy extrapolated x loc position at middle
955  sigma = abs(res) / sqrt(parameters.second.xx() + xErr * xErr / uxlden / uxlden +
956  yErr * yErr * xloc * xloc * uylfac * uylfac / uxlden / uxlden / uxlden / uxlden);
957  }
958 
959  VCluster_info.emplace_back(res, sigma, parameters.first.x());
960 
961  LogDebug("SiStripHitEfficiencyWorker") << "Have ID match. residual = " << res << " res sigma = " << sigma;
962  //LogDebug("SiStripHitEfficiencyWorker")
963  // << "trajectory measurement compatability estimate = " << (*itm).estimate() ;
964  LogDebug("SiStripHitEfficiencyWorker")
965  << "hit position = " << parameters.first.x() << " hit error = " << sqrt(parameters.second.xx())
966  << " trajectory position = " << xloc << " traj error = " << xErr;
967  }
968  }
969  ::ClusterInfo finalCluster{1000.0, 1000.0, 0.0};
970  if (!VCluster_info.empty()) {
971  LogDebug("SiStripHitEfficiencyWorker") << "found clusters > 0";
972  if (VCluster_info.size() > 1) {
973  //get the smallest one
974  for (const auto& res : VCluster_info) {
975  if (std::abs(res.xResidualPull) < std::abs(finalCluster.xResidualPull)) {
976  finalCluster = res;
977  }
978  LogDebug("SiStripHitEfficiencyWorker")
979  << "iresidual = " << res.xResidual << " isigma = " << res.xResidualPull
980  << " and FinalRes = " << finalCluster.xResidual;
981  }
982  } else {
983  finalCluster = VCluster_info[0];
984  }
985  VCluster_info.clear();
986  }
987 
988  LogDebug("SiStripHitEfficiencyWorker") << "Final residual in X = " << finalCluster.xResidual << "+-"
989  << (finalCluster.xResidual / finalCluster.xResidualPull);
990  LogDebug("SiStripHitEfficiencyWorker")
991  << "Checking location of trajectory: abs(yloc) = " << abs(yloc) << " abs(xloc) = " << abs(xloc);
992 
993  //
994  // fill ntuple varibles
995 
996  //if ( stripQuality->IsModuleBad(iidd) )
997  if (stripQuality.getBadApvs(iidd) != 0) {
998  SiStripQualBad = 1;
999  LogDebug("SiStripHitEfficiencyWorker") << "strip is bad from SiStripQuality";
1000  } else {
1001  SiStripQualBad = 0;
1002  LogDebug("SiStripHitEfficiencyWorker") << "strip is good from SiStripQuality";
1003  }
1004 
1005  //check for FED-detected errors and include those in SiStripQualBad
1006  for (unsigned int ii = 0; ii < fedErrorIds.size(); ii++) {
1007  if (iidd == fedErrorIds[ii].rawId())
1008  SiStripQualBad = 1;
1009  }
1010 
1011  // CM of APV crossed by traj
1012  if (addCommonMode_)
1013  if (commonModeDigis.isValid() && TrajStrip >= 0 && TrajStrip <= 768) {
1014  const auto digiframe = commonModeDigis->find(iidd);
1015  if (digiframe != commonModeDigis->end())
1016  if ((unsigned)TrajStrip / sistrip::STRIPS_PER_APV < digiframe->data.size())
1017  commonMode = digiframe->data.at(TrajStrip / sistrip::STRIPS_PER_APV).adc();
1018  }
1019 
1020  LogDebug("SiStripHitEfficiencyWorker") << "before check good";
1021 
1022  if (finalCluster.xResidualPull < 999.0) { //could make requirement on track/hit consistency, but for
1023  //now take anything with a hit on the module
1024  LogDebug("SiStripHitEfficiencyWorker")
1025  << "hit being counted as good " << finalCluster.xResidual << " FinalRecHit " << iidd << " TKlayers "
1026  << TKlayers << " xloc " << xloc << " yloc " << yloc << " module " << iidd
1027  << " matched/stereo/rphi = " << ((iidd & 0x3) == 0) << "/" << ((iidd & 0x3) == 1) << "/"
1028  << ((iidd & 0x3) == 2);
1029  ModIsBad = 0;
1030  } else {
1031  LogDebug("SiStripHitEfficiencyWorker")
1032  << "hit being counted as bad ######### Invalid RPhi FinalResX " << finalCluster.xResidual
1033  << " FinalRecHit " << iidd << " TKlayers " << TKlayers << " xloc " << xloc << " yloc " << yloc
1034  << " module " << iidd << " matched/stereo/rphi = " << ((iidd & 0x3) == 0) << "/" << ((iidd & 0x3) == 1)
1035  << "/" << ((iidd & 0x3) == 2);
1036  ModIsBad = 1;
1037  LogDebug("SiStripHitEfficiencyWorker")
1038  << " RPhi Error " << sqrt(xErr * xErr + yErr * yErr) << " ErrorX " << xErr << " yErr " << yErr;
1039  }
1040 
1041  LogDebug("SiStripHitEfficiencyWorker")
1042  << "To avoid them staying unused: ModIsBad=" << ModIsBad << ", SiStripQualBad=" << SiStripQualBad
1043  << ", commonMode=" << commonMode << ", highPurity=" << highPurity
1044  << ", withinAcceptance=" << withinAcceptance;
1045 
1046  unsigned int layer = TKlayers;
1047  if (showRings_ && layer > bounds::k_LayersAtTOBEnd) { // use rings instead of wheels
1048  if (layer <= bounds::k_LayersAtTIDEnd) { // TID
1049  layer = bounds::k_LayersAtTOBEnd +
1050  tTopo->tidRing(iidd); // ((iidd >> 9) & 0x3); // 3 disks and also 3 rings -> use the same container
1051  } else { // TEC
1052  layer = bounds::k_LayersAtTIDEnd + tTopo->tecRing(iidd); // ((iidd >> 5) & 0x7);
1053  }
1054  }
1055  unsigned int layerWithSide = layer;
1056  if (layer > bounds::k_LayersAtTOBEnd && layer <= bounds::k_LayersAtTIDEnd) {
1057  const auto side = tTopo->tidSide(iidd); //(iidd >> 13) & 0x3; // TID
1058  if (side == 2)
1059  layerWithSide = layer + 3;
1060  } else if (layer > bounds::k_LayersAtTIDEnd) {
1061  const auto side = tTopo->tecSide(iidd); // (iidd >> 18) & 0x3; // TEC
1062  if (side == 1) {
1063  layerWithSide = layer + 3;
1064  } else if (side == 2) {
1065  layerWithSide = layer + 3 + (showRings_ ? 7 : 9);
1066  }
1067  }
1068 
1069  if ((bunchX_ > 0 && bunchX_ != bunchCrossing) || (!withinAcceptance) ||
1071  (!showTOB6TEC9_ && (TKlayers == bounds::k_LayersAtTOBEnd || TKlayers == bounds::k_LayersAtTECEnd)) ||
1072  (badModules_.end() != badModules_.find(iidd)))
1073  return;
1074 
1075  const bool badquality = (SiStripQualBad == 1);
1076 
1077  //Now that we have a good event, we need to look at if we expected it or not, and the location
1078  //if we didn't
1079  //Fill the missing hit information first
1080  bool badflag = false; // true for hits that are expected but not found
1081  if (resXSig_ < 0) {
1082  if (ModIsBad == 1)
1083  badflag = true; // isBad set to false in the tree when resxsig<999.0
1084  } else {
1085  if (ModIsBad == 1 || finalCluster.xResidualPull > resXSig_)
1086  badflag = true;
1087  }
1088 
1089  // Conversion of positions in strip unit
1090  int nstrips = -9;
1091  float Pitch = -9.0;
1092  const StripGeomDetUnit* stripdet = nullptr;
1093  if (finalCluster.xResidualPull ==
1094  1000.0) { // special treatment, no GeomDetUnit associated in some cases when no cluster found
1095  Pitch = 0.0205; // maximum
1096  nstrips = 768; // maximum
1097  } else {
1098  stripdet = dynamic_cast<const StripGeomDetUnit*>(tkgeom->idToDetUnit(iidd));
1099  const StripTopology& Topo = stripdet->specificTopology();
1100  nstrips = Topo.nstrips();
1101  Pitch = stripdet->surface().bounds().width() / Topo.nstrips();
1102  }
1103  double stripTrajMid = xloc / Pitch + nstrips / 2.0;
1104  double stripCluster = finalCluster.xLocal / Pitch + nstrips / 2.0;
1105  // For trapezoidal modules: extrapolation of x trajectory position to the y middle of the module
1106  // for correct comparison with cluster position
1107  if (stripdet && layer > bounds::k_LayersAtTOBEnd) {
1108  const auto& trapezoidalBounds = dynamic_cast<const TrapezoidalPlaneBounds&>(stripdet->surface().bounds());
1109  std::array<const float, 4> const& parameters = trapezoidalBounds.parameters();
1110  const float hbedge = parameters[0];
1111  const float htedge = parameters[1];
1112  const float hapoth = parameters[3];
1113  const float TrajLocXMid = xloc / (1 + (htedge - hbedge) * yloc / (htedge + hbedge) /
1114  hapoth); // radialy extrapolated x loc position at middle
1115  stripTrajMid = TrajLocXMid / Pitch + nstrips / 2.0;
1116  }
1117 
1118  if ((!badquality) && (layer < h_resolution.size())) {
1119  LogDebug("SiStripHitEfficiencyWorker")
1120  << "layer " << layer << " vector index " << layer - 1 << " before filling h_resolution" << std::endl;
1121  h_resolution[layer - 1]->Fill(finalCluster.xResidualPull != 1000.0 ? stripTrajMid - stripCluster : 1000);
1122  }
1123 
1124  // New matching methods
1125  if (clusterMatchingMethod_ >= 1) {
1126  badflag = false;
1127  if (finalCluster.xResidualPull == 1000.0) {
1128  LogDebug("SiStripHitEfficiencyWorker") << "Marking bad for resxsig=1000";
1129  badflag = true;
1130  } else {
1132  // check the distance between cluster and trajectory position
1133  if (std::abs(stripCluster - stripTrajMid) > clusterTracjDist_) {
1134  LogDebug("SiStripHitEfficiencyWorker") << "Marking bad for cluster-to-traj distance";
1135  badflag = true;
1136  }
1137  }
1139  // cluster and traj have to be in the same APV (don't take edges into accounts)
1140  const int tapv = (int)stripTrajMid / sistrip::STRIPS_PER_APV;
1141  const int capv = (int)stripCluster / sistrip::STRIPS_PER_APV;
1142  float stripInAPV = stripTrajMid - tapv * sistrip::STRIPS_PER_APV;
1143  if (stripInAPV < stripsApvEdge_ || stripInAPV > sistrip::STRIPS_PER_APV - stripsApvEdge_) {
1144  LogDebug("SiStripHitEfficiencyWorker") << "Too close to the edge: " << stripInAPV;
1145  return;
1146  }
1147  if (tapv != capv) {
1148  LogDebug("SiStripHitEfficiencyWorker") << "Marking bad for tapv!=capv";
1149  badflag = true;
1150  }
1151  }
1152  }
1153  }
1154  if (!badquality) {
1155  LogDebug("SiStripHitEfficiencyWorker")
1156  << "Filling measurement for " << iidd << " in layer " << layer << " histograms with bx=" << bunchCrossing
1157  << ", lumi=" << instLumi << ", PU=" << PU << "; bad flag=" << badflag;
1158 
1159  // hot/cold maps of hits that are expected but not found
1160  if (badflag) {
1161  if (layer > bounds::k_LayersStart && layer <= bounds::k_LayersAtTIBEnd) {
1162  //We are in the TIB
1163  float phi = ::calcPhi(tm.globalX(), tm.globalY());
1164  h_hotcold[layer - 1]->Fill(360. - phi, tm.globalZ(), 1.);
1165  } else if (layer > bounds::k_LayersAtTIBEnd && layer <= bounds::k_LayersAtTOBEnd) {
1166  //We are in the TOB
1167  float phi = ::calcPhi(tm.globalX(), tm.globalY());
1168  h_hotcold[layer - 1]->Fill(360. - phi, tm.globalZ(), 1.);
1169  } else if (layer > bounds::k_LayersAtTOBEnd && layer <= bounds::k_LayersAtTIDEnd) {
1170  //We are in the TID
1171  //There are 2 different maps here
1172  int side = tTopo->tidSide(iidd);
1173  if (side == 1)
1174  h_hotcold[(layer - 1) + (layer - 11)]->Fill(-tm.globalY(), tm.globalX(), 1.);
1175  else if (side == 2)
1176  h_hotcold[(layer - 1) + (layer - 10)]->Fill(-tm.globalY(), tm.globalX(), 1.);
1177  } else if (layer > bounds::k_LayersAtTIDEnd) {
1178  //We are in the TEC
1179  //There are 2 different maps here
1180  int side = tTopo->tecSide(iidd);
1181  if (side == 1)
1182  h_hotcold[(layer + 2) + (layer - 14)]->Fill(-tm.globalY(), tm.globalX(), 1.);
1183  else if (side == 2)
1184  h_hotcold[(layer + 2) + (layer - 13)]->Fill(-tm.globalY(), tm.globalX(), 1.);
1185  }
1186  }
1187 
1188  LogDebug("SiStripHitEfficiencyWorker")
1189  << "layer " << layer << " vector index " << layer - 1 << " before filling h_layer_vsSmthg" << std::endl;
1190  h_layer_vsBx[layer - 1].fill(bunchCrossing, !badflag);
1191  if (addLumi_) {
1192  h_layer_vsLumi[layer - 1].fill(instLumi, !badflag);
1193  h_layer_vsPU[layer - 1].fill(PU, !badflag);
1194  }
1195  if (addCommonMode_) {
1196  h_layer_vsCM[layer - 1].fill(commonMode, !badflag);
1197  }
1198  h_goodLayer.fill(layerWithSide, !badflag);
1199  }
1200  // efficiency without bad modules excluded
1201  h_allLayer.fill(layerWithSide, !badflag);
1202 
1203  // efficiency without bad modules excluded
1204  if (TKlayers) {
1205  h_module.fill(iidd, !badflag);
1206  assert(h_module.check(iidd));
1207  }
1208 
1209  /* Used in SiStripHitEffFromCalibTree:
1210  * run -> "run" -> run // e.id().run()
1211  * event -> "event" -> evt // e.id().event()
1212  * ModIsBad -> "ModIsBad" -> isBad
1213  * SiStripQualBad -> "SiStripQualBad"" -> quality
1214  * Id -> "Id" -> id // iidd
1215  * withinAcceptance -> "withinAcceptance" -> accept
1216  * whatlayer -> "layer" -> layer_wheel // Tklayers
1217  * highPurity -> "highPurity" -> highPurity
1218  * TrajGlbX -> "TrajGlbX" -> x // tm.globalX()
1219  * TrajGlbY -> "TrajGlbY" -> y // tm.globalY()
1220  * TrajGlbZ -> "TrajGlbZ" -> z // tm.globalZ()
1221  * ResXSig -> "ResXSig" -> resxsig // finalCluster.xResidualPull;
1222  * TrajLocX -> "TrajLocX" -> TrajLocX // xloc
1223  * TrajLocY -> "TrajLocY" -> TrajLocY // yloc
1224  * ClusterLocX -> "ClusterLocX" -> ClusterLocX // finalCluster.xLocal
1225  * bunchx -> "bunchx" -> bx // e.bunchCrossing()
1226  * instLumi -> "instLumi" -> instLumi ## if addLumi_
1227  * PU -> "PU" -> PU ## if addLumi_
1228  * commonMode -> "commonMode" -> CM ## if addCommonMode_ / _useCM
1229  */
1230  LogDebug("SiStripHitEfficiencyWorker") << "after good location check";
1231  }
1232  LogDebug("SiStripHitEfficiencyWorker") << "after list of clusters";
1233  }
1234  LogDebug("SiStripHitEfficiencyWorker") << "After layers=TKLayers if with TKlayers=" << TKlayers
1235  << ", layers=" << layers_;
1236 }
1237 
1240  desc.add<std::string>("dqmDir", "AlCaReco/SiStripHitEfficiency");
1241  desc.add<bool>("UseOnlyHighPurityTracks", true);
1242  desc.add<bool>("cutOnTracks", false);
1243  desc.add<bool>("doMissingHitsRecovery", false);
1244  desc.add<bool>("useAllHitsFromTracksWithMissingHits", false);
1245  desc.add<bool>("useFirstMeas", false);
1246  desc.add<bool>("useLastMeas", false);
1247  desc.add<double>("ClusterTrajDist", 64.0);
1248  desc.add<double>("ResXSig", -1);
1249  desc.add<double>("StripsApvEdge", 10.0);
1250  desc.add<edm::InputTag>("combinatorialTracks", edm::InputTag{"generalTracks"});
1251  desc.add<edm::InputTag>("commonMode", edm::InputTag{"siStripDigis", "CommonMode"});
1252  desc.add<edm::InputTag>("lumiScalers", edm::InputTag{"scalersRawToDigi"});
1253  desc.add<edm::InputTag>("metadata", edm::InputTag{"onlineMetaDataDigis"});
1254  desc.add<edm::InputTag>("siStripClusters", edm::InputTag{"siStripClusters"});
1255  desc.add<edm::InputTag>("siStripDigis", edm::InputTag{"siStripDigis"});
1256  desc.add<edm::InputTag>("trackerEvent", edm::InputTag{"MeasurementTrackerEvent"});
1257  desc.add<edm::InputTag>("trajectories", edm::InputTag{"generalTracks"});
1258  desc.add<int>("ClusterMatchingMethod", 0);
1259  desc.add<int>("Layer", 0);
1260  desc.add<unsigned int>("trackMultiplicity", 100);
1261  desc.addUntracked<bool>("Debug", false);
1262  desc.addUntracked<bool>("ShowRings", false);
1263  desc.addUntracked<bool>("ShowTOB6TEC9", false);
1264  desc.addUntracked<bool>("addCommonMode", false);
1265  desc.addUntracked<bool>("addLumi", true);
1266  desc.addUntracked<int>("BunchCrossing", 0);
1267  desc.addUntracked<std::string>("BadModulesFile", "");
1268  descriptions.addWithDefaultLabel(desc);
1269 }
1270 
std::pair< LocalPoint, LocalError > LocalValues
const edm::ESGetToken< MeasurementTracker, CkfComponentsRecord > measTrackerToken_
const edm::EDGetTokenT< DetIdVector > digisVec_token_
virtual int nstrips() const =0
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
void fill(double x, bool found, float weight=1.)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
virtual std::vector< DetWithState > compatibleDets(const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est) const
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
std::vector< MonitorElement * > h_resolution
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
EffME1(MonitorElement *total, MonitorElement *found)
void fillForTraj(const TrajectoryAtInvalidHit &tm, const TrackerTopology *tTopo, const TrackerGeometry *tkgeom, const StripClusterParameterEstimator &stripCPE, const SiStripQuality &stripQuality, const DetIdVector &fedErrorIds, const edm::Handle< edm::DetSetVector< SiStripRawDigi >> &commonModeDigis, const edmNew::DetSetVector< SiStripCluster > &theClusters, int bunchCrossing, float instLumi, float PU, bool highPurity)
const edm::EDGetTokenT< std::vector< Trajectory > > trajectories_token_
std::vector< MonitorElement * > h_hotcold
virtual void localParameters(AClusters const &clusters, ALocalValues &retValues, const GeomDetUnit &gd, const LocalTrajectoryParameters &ltp) const
void bookHistograms(DQMStore::IBooker &booker, const edm::Run &run, const edm::EventSetup &setup) override
unsigned int tidSide(const DetId &id) const
Class to contain the online luminosity from soft FED 1022.
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
unsigned int monodet_id() const
Definition: weight.py:1
float instLumi() const
Return the luminosity for the current nibble.
std::unordered_map< uint32_t, int > fedErrorCounts
assert(be >=bs)
unsigned int tecRing(const DetId &id) const
ring id
Definition: Electron.h:6
missing
Definition: combine.py:5
std::unique_ptr< TkHistoMap > FEDErrorOccupancy
const edm::EDGetTokenT< MeasurementTrackerEvent > trackerEvent_token_
const_iterator end(bool update=false) const
~SiStripHitEfficiencyWorker() override=default
void Fill(long long x)
const edm::ESGetToken< TkDetMap, TrackerTopologyRcd > tkDetMapToken_
const edm::ESGetToken< StripClusterParameterEstimator, TkStripCPERecord > stripCPEToken_
const edm::ESGetToken< Propagator, TrackingComponentsRecord > propagatorToken_
const edm::EDGetTokenT< OnlineLuminosityRecord > metaDataToken_
const edm::ESGetToken< Chi2MeasurementEstimatorBase, TrackingComponentsRecord > chi2EstimatorToken_
unsigned int tecSide(const DetId &id) const
EffTkMap(std::unique_ptr< TkHistoMap > &&total, std::unique_ptr< TkHistoMap > &&found)
MonitorElement * bookProfile(TString const &name, TString const &title, int nchX, double lowX, double highX, int, double lowY, double highY, char const *option="s", FUNC onbooking=NOOP())
Definition: DQMStore.h:399
T sqrt(T t)
Definition: SSEVec.h:19
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
virtual const StripTopology & specificTopology() const
Returns a reference to the strip proxy topology.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
virtual void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
dqm::reco::MonitorElement * EventStats
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
ii
Definition: cuy.py:589
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magFieldToken_
SiStripHitEfficiencyWorker(const edm::ParameterSet &conf)
void analyze(const edm::Event &e, const edm::EventSetup &c) override
std::vector< DetId > DetIdVector
Definition: DetIdVector.h:7
std::vector< unsigned int > hitTotalCounters
Definition: DetId.h:17
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
size_type size() const
const edm::EDGetTokenT< TrajTrackAssociationCollection > trajTrackAsso_token_
Constants and enumerated types for FED/FEC systems.
part
Definition: HCALResponse.h:20
virtual const std::array< const float, 4 > parameters() const
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:212
const edm::EDGetTokenT< DetIdCollection > digisCol_token_
const edm::EDGetTokenT< edmNew::DetSetVector< SiStripCluster > > clusters_token_
bool isValid() const
Definition: HandleBase.h:70
fixed size matrix
HLT enums.
const_iterator find(id_type i, bool update=false) const
static const uint16_t STRIPS_PER_APV
const edm::EDGetTokenT< reco::TrackCollection > combinatorialTracks_token_
const edm::ESGetToken< SiStripQuality, SiStripQualityRcd > stripQualityToken_
const edm::EDGetTokenT< LumiScalersCollection > scalerToken_
std::vector< unsigned int > hitRecoveryCounters
unsigned int tidRing(const DetId &id) const
std::vector< LumiScalers > LumiScalersCollection
Definition: LumiScalers.h:144
void fill(uint32_t id, bool found, float weight=1.)
const edm::EDGetTokenT< edm::DetSetVector< SiStripRawDigi > > commonModeToken_
Log< level::Warning, false > LogWarning
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
A Digi for the silicon strip detector, containing only adc information, and suitable for storing raw ...
tmp
align.sh
Definition: createJobs.py:716
virtual float width() const =0
MonitorElement * book2I(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:296
float avgPileUp() const
Return the average pileup for th current nibble.
def move(src, dest)
Definition: eostools.py:511
Definition: Run.h:45
ConstRecHitPointer const & recHit() const
#define LogDebug(id)
virtual void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
short getBadApvs(uint32_t detid) const
const Bounds & bounds() const
Definition: Surface.h:87