CMS 3D CMS Logo

MtdTracksValidation.cc
Go to the documentation of this file.
1 #define EDM_ML_DEBUG
2 #include <string>
3 
9 
12 
19 
25 
29 
37 
48 
59 
62 
69 
72 #include "HepMC/GenRanges.h"
73 #include "CLHEP/Units/PhysicalConstants.h"
74 #include "MTDHit.h"
75 
77 public:
78  explicit MtdTracksValidation(const edm::ParameterSet&);
79  ~MtdTracksValidation() override;
80 
81  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
82 
83 private:
84  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
85 
86  void analyze(const edm::Event&, const edm::EventSetup&) override;
87 
88  const std::pair<bool, bool> checkAcceptance(
89  const reco::Track&, const edm::Event&, const edm::EventSetup&, size_t&, float&, float&, float&, float&);
90 
91  const bool mvaGenSel(const HepMC::GenParticle&, const float&);
92  const bool mvaTPSel(const TrackingParticle&);
93  const bool mvaRecSel(const reco::TrackBase&, const reco::Vertex&, const double&, const double&);
94  const bool mvaGenRecMatch(const HepMC::GenParticle&, const double&, const reco::TrackBase&, const bool&);
96 
97  const unsigned long int uniqueId(const uint32_t x, const EncodedEventId& y) {
98  const uint64_t a = static_cast<uint64_t>(x);
99  const uint64_t b = static_cast<uint64_t>(y.rawId());
100 
101  if (x < y.rawId())
102  return (b << 32) | a;
103  else
104  return (a << 32) | b;
105  }
106 
107  bool isETL(const double eta) const { return (std::abs(eta) > trackMinEtlEta_) && (std::abs(eta) < trackMaxEtlEta_); }
108 
109  // ------------ member data ------------
110 
112  const float trackMinPt_;
113  const float trackMaxBtlEta_;
114  const float trackMinEtlEta_;
115  const float trackMaxEtlEta_;
116 
117  static constexpr double etacutGEN_ = 4.; // |eta| < 4;
118  static constexpr double etacutREC_ = 3.; // |eta| < 3;
119  static constexpr double pTcut_ = 0.7; // PT > 0.7 GeV
120  static constexpr double deltaZcut_ = 0.1; // dz separation 1 mm
121  static constexpr double deltaPTcut_ = 0.05; // dPT < 5%
122  static constexpr double deltaDRcut_ = 0.03; // DeltaR separation
123  static constexpr double depositBTLthreshold_ = 1; // threshold for energy deposit in BTL cell [MeV]
124  static constexpr double depositETLthreshold_ = 0.001; // threshold for energy deposit in ETL cell [MeV]
125  static constexpr double rBTL_ = 110.0;
126  static constexpr double zETL_ = 290.0;
127  static constexpr double etaMatchCut_ = 0.05;
128  static constexpr double cluDRradius_ = 0.05; // to cluster rechits around extrapolated track
129 
132 
136 
145 
148 
158 
165 
174 
186 
198 
205 
224 
244 
250 };
251 
252 // ------------ constructor and destructor --------------
254  : folder_(iConfig.getParameter<std::string>("folder")),
255  trackMinPt_(iConfig.getParameter<double>("trackMinimumPt")),
256  trackMaxBtlEta_(iConfig.getParameter<double>("trackMaximumBtlEta")),
257  trackMinEtlEta_(iConfig.getParameter<double>("trackMinimumEtlEta")),
258  trackMaxEtlEta_(iConfig.getParameter<double>("trackMaximumEtlEta")) {
259  GenRecTrackToken_ = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("inputTagG"));
260  RecTrackToken_ = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("inputTagT"));
261  RecVertexToken_ = consumes<std::vector<reco::Vertex>>(iConfig.getParameter<edm::InputTag>("inputTagV"));
262  HepMCProductToken_ = consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("inputTagH"));
264  consumes<TrackingParticleCollection>(iConfig.getParameter<edm::InputTag>("SimTag"));
266  consumes<reco::SimToRecoCollection>(iConfig.getParameter<edm::InputTag>("TPtoRecoTrackAssoc"));
268  consumes<reco::RecoToSimCollection>(iConfig.getParameter<edm::InputTag>("TPtoRecoTrackAssoc"));
269  btlSimHitsToken_ = consumes<CrossingFrame<PSimHit>>(iConfig.getParameter<edm::InputTag>("btlSimHits"));
270  etlSimHitsToken_ = consumes<CrossingFrame<PSimHit>>(iConfig.getParameter<edm::InputTag>("etlSimHits"));
271  btlRecHitsToken_ = consumes<FTLRecHitCollection>(iConfig.getParameter<edm::InputTag>("btlRecHits"));
272  etlRecHitsToken_ = consumes<FTLRecHitCollection>(iConfig.getParameter<edm::InputTag>("etlRecHits"));
273  trackAssocToken_ = consumes<edm::ValueMap<int>>(iConfig.getParameter<edm::InputTag>("trackAssocSrc"));
274  pathLengthToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("pathLengthSrc"));
275  tmtdToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("tmtd"));
276  SigmatmtdToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmatmtd"));
277  t0SrcToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("t0Src"));
278  Sigmat0SrcToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmat0Src"));
279  t0PidToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("t0PID"));
280  Sigmat0PidToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmat0PID"));
281  t0SafePidToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("t0SafePID"));
282  Sigmat0SafePidToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmat0SafePID"));
283  trackMVAQualToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("trackMVAQual"));
284  mtdgeoToken_ = esConsumes<MTDGeometry, MTDDigiGeometryRecord>();
285  mtdtopoToken_ = esConsumes<MTDTopology, MTDTopologyRcd>();
286  mtdlayerToken_ = esConsumes<MTDDetLayerGeometry, MTDRecoGeometryRecord>();
287  magfieldToken_ = esConsumes<MagneticField, IdealMagneticFieldRecord>();
288  builderToken_ = esConsumes<TransientTrackBuilder, TransientTrackRecord>(edm::ESInputTag("", "TransientTrackBuilder"));
289  particleTableToken_ = esConsumes<HepPDT::ParticleDataTable, edm::DefaultRecord>();
290 }
291 
293 
294 // ------------ method called for each event ------------
296  using namespace edm;
297  using namespace geant_units::operators;
298  using namespace std;
299 
300  auto GenRecTrackHandle = makeValid(iEvent.getHandle(GenRecTrackToken_));
301  auto RecVertexHandle = makeValid(iEvent.getHandle(RecVertexToken_));
302 
303  std::unordered_map<uint32_t, MTDHit> m_btlHits;
304  std::unordered_map<uint32_t, MTDHit> m_etlHits;
305  std::unordered_map<uint32_t, std::set<unsigned long int>> m_btlTrkPerCell;
306  std::unordered_map<uint32_t, std::set<unsigned long int>> m_etlTrkPerCell;
307  std::map<TrackingParticleRef, std::vector<uint32_t>> m_tp2detid;
308 
309  const auto& tMtd = iEvent.get(tmtdToken_);
310  const auto& SigmatMtd = iEvent.get(SigmatmtdToken_);
311  const auto& t0Src = iEvent.get(t0SrcToken_);
312  const auto& Sigmat0Src = iEvent.get(Sigmat0SrcToken_);
313  const auto& t0Pid = iEvent.get(t0PidToken_);
314  const auto& Sigmat0Pid = iEvent.get(Sigmat0PidToken_);
315  const auto& t0Safe = iEvent.get(t0SafePidToken_);
316  const auto& Sigmat0Safe = iEvent.get(Sigmat0SafePidToken_);
317  const auto& mtdQualMVA = iEvent.get(trackMVAQualToken_);
318  const auto& trackAssoc = iEvent.get(trackAssocToken_);
319  const auto& pathLength = iEvent.get(pathLengthToken_);
320 
321  const auto& primRecoVtx = *(RecVertexHandle.product()->begin());
322 
323  // generator level information (HepMC format)
324  auto GenEventHandle = makeValid(iEvent.getHandle(HepMCProductToken_));
325  const HepMC::GenEvent* mc = GenEventHandle->GetEvent();
326  double zsim = convertMmToCm((*(mc->vertices_begin()))->position().z());
327  double tsim = (*(mc->vertices_begin()))->position().t() * CLHEP::mm / CLHEP::c_light;
328 
329  auto pdt = iSetup.getHandle(particleTableToken_);
330  const HepPDT::ParticleDataTable* pdTable = pdt.product();
331 
332  auto simToRecoH = makeValid(iEvent.getHandle(simToRecoAssociationToken_));
333  s2r_ = simToRecoH.product();
334 
335  auto recoToSimH = makeValid(iEvent.getHandle(recoToSimAssociationToken_));
336  r2s_ = recoToSimH.product();
337 
338  //Fill maps with simhits accumulated per DetId
339 
340  auto btlSimHitsHandle = makeValid(iEvent.getHandle(btlSimHitsToken_));
341  MixCollection<PSimHit> btlSimHits(btlSimHitsHandle.product());
342  for (auto const& simHit : btlSimHits) {
343  if (simHit.tof() < 0 || simHit.tof() > 25.)
344  continue;
345  DetId id = simHit.detUnitId();
346  auto const thisHId = uniqueId(simHit.trackId(), simHit.eventId());
347  m_btlTrkPerCell[id.rawId()].insert(thisHId);
348  auto simHitIt = m_btlHits.emplace(id.rawId(), MTDHit()).first;
349  // --- Accumulate the energy (in MeV) of SIM hits in the same detector cell
350  (simHitIt->second).energy += convertUnitsTo(0.001_MeV, simHit.energyLoss());
351  }
352 
353  auto etlSimHitsHandle = makeValid(iEvent.getHandle(etlSimHitsToken_));
354  MixCollection<PSimHit> etlSimHits(etlSimHitsHandle.product());
355  for (auto const& simHit : etlSimHits) {
356  if (simHit.tof() < 0 || simHit.tof() > 25.) {
357  continue;
358  }
359  DetId id = simHit.detUnitId();
360  auto const thisHId = uniqueId(simHit.trackId(), simHit.eventId());
361  m_etlTrkPerCell[id.rawId()].insert(thisHId);
362  auto simHitIt = m_etlHits.emplace(id.rawId(), MTDHit()).first;
363  // --- Accumulate the energy (in MeV) of SIM hits in the same detector cell
364  (simHitIt->second).energy += convertUnitsTo(0.001_MeV, simHit.energyLoss());
365  }
366 
367  //Fill map of DetId per ref to TP
368 
369  auto tpHandle = makeValid(iEvent.getHandle(trackingParticleCollectionToken_));
370  TrackingParticleCollection tpColl = *(tpHandle.product());
371  size_t tpindex(0);
372  for (auto tp = tpColl.begin(); tp != tpColl.end(); tp++, ++tpindex) {
374  if (tp->eventId().bunchCrossing() == 0 && tp->eventId().event() == 0) {
375  if (!mvaTPSel(*tp))
376  continue;
377  for (const auto& simTrk : tp->g4Tracks()) {
378  auto const thisTId = uniqueId(simTrk.trackId(), simTrk.eventId());
379  for (auto const& cell : m_btlTrkPerCell) {
380  if (m_btlHits[cell.first].energy < depositBTLthreshold_) {
381  continue;
382  }
383  for (auto const& simtrack : cell.second) {
384  if (thisTId == simtrack) {
385  m_tp2detid[tpref].emplace_back(cell.first);
386  break;
387  }
388  }
389  }
390  for (auto const& cell : m_etlTrkPerCell) {
391  if (m_etlHits[cell.first].energy < depositETLthreshold_) {
392  continue;
393  }
394  for (auto const& simtrack : cell.second) {
395  if (thisTId == simtrack) {
396  m_tp2detid[tpref].emplace_back(cell.first);
397  break;
398  }
399  }
400  }
401  }
402  }
403  }
404 
405  unsigned int index = 0;
406 
407  // flag to select events with reco vertex close to true simulated primary vertex, or PV fake (particle guns)
408  const bool isGoodVtx = std::abs(primRecoVtx.z() - zsim) < deltaZcut_ || primRecoVtx.isFake();
409 
410  // --- Loop over all RECO tracks ---
411  for (const auto& trackGen : *GenRecTrackHandle) {
412  const reco::TrackRef trackref(iEvent.getHandle(GenRecTrackToken_), index);
413  index++;
414 
415  if (trackAssoc[trackref] == -1) {
416  LogInfo("mtdTracks") << "Extended track not associated";
417  continue;
418  }
419 
420  const reco::TrackRef mtdTrackref = reco::TrackRef(iEvent.getHandle(RecTrackToken_), trackAssoc[trackref]);
421  const reco::Track& track = *mtdTrackref;
422 
423  bool isBTL = false;
424  bool isETL = false;
425  bool twoETLdiscs = false;
426  bool noCrack = std::abs(trackGen.eta()) < trackMaxBtlEta_ || std::abs(trackGen.eta()) > trackMinEtlEta_;
427 
428  if (track.pt() >= trackMinPt_ && std::abs(track.eta()) <= trackMaxEtlEta_) {
429  meTracktmtd_->Fill(tMtd[trackref]);
430  if (std::round(SigmatMtd[trackref] - Sigmat0Pid[trackref]) != 0) {
431  LogWarning("mtdTracks")
432  << "TimeError associated to refitted track is different from TimeError stored in tofPID "
433  "sigmat0 ValueMap: this should not happen";
434  }
435 
436  meTrackt0Src_->Fill(t0Src[trackref]);
437  meTrackSigmat0Src_->Fill(Sigmat0Src[trackref]);
438 
439  meTrackt0Pid_->Fill(t0Pid[trackref]);
440  meTrackSigmat0Pid_->Fill(Sigmat0Pid[trackref]);
441  meTrackt0SafePid_->Fill(t0Safe[trackref]);
442  meTrackSigmat0SafePid_->Fill(Sigmat0Safe[trackref]);
443  meTrackMVAQual_->Fill(mtdQualMVA[trackref]);
444 
445  meTrackPathLenghtvsEta_->Fill(std::abs(track.eta()), pathLength[trackref]);
446 
447  if (std::abs(track.eta()) < trackMaxBtlEta_) {
448  // --- all BTL tracks (with and without hit in MTD) ---
452 
453  bool MTDBtl = false;
454  int numMTDBtlvalidhits = 0;
455  for (const auto hit : track.recHits()) {
456  if (hit->isValid() == false)
457  continue;
458  MTDDetId Hit = hit->geographicalId();
459  if ((Hit.det() == 6) && (Hit.subdetId() == 1) && (Hit.mtdSubDetector() == 1)) {
460  MTDBtl = true;
461  numMTDBtlvalidhits++;
462  }
463  }
464  meTrackNumHits_->Fill(numMTDBtlvalidhits);
465 
466  // --- keeping only tracks with last hit in MTD ---
467  if (MTDBtl == true) {
468  isBTL = true;
473  meBTLTrackPtRes_->Fill((trackGen.pt() - track.pt()) / trackGen.pt());
474  }
475  if (isBTL && Sigmat0Safe[trackref] < 0.) {
476  meTrackNumHitsNT_->Fill(numMTDBtlvalidhits);
477  }
478  } //loop over (geometrical) BTL tracks
479 
480  else {
481  // --- all ETL tracks (with and without hit in MTD) ---
482  if ((track.eta() < -trackMinEtlEta_) && (track.eta() > -trackMaxEtlEta_)) {
483  meETLTrackEffEtaTot_[0]->Fill(track.eta());
484  meETLTrackEffPhiTot_[0]->Fill(track.phi());
485  meETLTrackEffPtTot_[0]->Fill(track.pt());
486  }
487 
488  if ((track.eta() > trackMinEtlEta_) && (track.eta() < trackMaxEtlEta_)) {
489  meETLTrackEffEtaTot_[1]->Fill(track.eta());
490  meETLTrackEffPhiTot_[1]->Fill(track.phi());
491  meETLTrackEffPtTot_[1]->Fill(track.pt());
492  }
493 
494  bool MTDEtlZnegD1 = false;
495  bool MTDEtlZnegD2 = false;
496  bool MTDEtlZposD1 = false;
497  bool MTDEtlZposD2 = false;
498  int numMTDEtlvalidhits = 0;
499  for (const auto hit : track.recHits()) {
500  if (hit->isValid() == false)
501  continue;
502  MTDDetId Hit = hit->geographicalId();
503  if ((Hit.det() == 6) && (Hit.subdetId() == 1) && (Hit.mtdSubDetector() == 2)) {
504  isETL = true;
505  ETLDetId ETLHit = hit->geographicalId();
506 
507  if ((ETLHit.zside() == -1) && (ETLHit.nDisc() == 1)) {
508  MTDEtlZnegD1 = true;
510  meETLTrackPtRes_->Fill((trackGen.pt() - track.pt()) / trackGen.pt());
511  numMTDEtlvalidhits++;
512  }
513  if ((ETLHit.zside() == -1) && (ETLHit.nDisc() == 2)) {
514  MTDEtlZnegD2 = true;
516  meETLTrackPtRes_->Fill((trackGen.pt() - track.pt()) / trackGen.pt());
517  numMTDEtlvalidhits++;
518  }
519  if ((ETLHit.zside() == 1) && (ETLHit.nDisc() == 1)) {
520  MTDEtlZposD1 = true;
522  meETLTrackPtRes_->Fill((trackGen.pt() - track.pt()) / trackGen.pt());
523  numMTDEtlvalidhits++;
524  }
525  if ((ETLHit.zside() == 1) && (ETLHit.nDisc() == 2)) {
526  MTDEtlZposD2 = true;
528  meETLTrackPtRes_->Fill((trackGen.pt() - track.pt()) / trackGen.pt());
529  numMTDEtlvalidhits++;
530  }
531  }
532  }
533  meTrackNumHits_->Fill(-numMTDEtlvalidhits);
534  if (isETL && Sigmat0Safe[trackref] < 0.) {
535  meTrackNumHitsNT_->Fill(-numMTDEtlvalidhits);
536  }
537 
538  // --- keeping only tracks with last hit in MTD ---
539  if ((track.eta() < -trackMinEtlEta_) && (track.eta() > -trackMaxEtlEta_)) {
540  twoETLdiscs = (MTDEtlZnegD1 == true) && (MTDEtlZnegD2 == true);
541  if ((MTDEtlZnegD1 == true) || (MTDEtlZnegD2 == true)) {
542  meETLTrackEffEtaMtd_[0]->Fill(track.eta());
543  meETLTrackEffPhiMtd_[0]->Fill(track.phi());
544  meETLTrackEffPtMtd_[0]->Fill(track.pt());
545  if (twoETLdiscs) {
546  meETLTrackEffEta2Mtd_[0]->Fill(track.eta());
547  meETLTrackEffPhi2Mtd_[0]->Fill(track.phi());
548  meETLTrackEffPt2Mtd_[0]->Fill(track.pt());
549  }
550  }
551  }
552  if ((track.eta() > trackMinEtlEta_) && (track.eta() < trackMaxEtlEta_)) {
553  twoETLdiscs = (MTDEtlZposD1 == true) && (MTDEtlZposD2 == true);
554  if ((MTDEtlZposD1 == true) || (MTDEtlZposD2 == true)) {
555  meETLTrackEffEtaMtd_[1]->Fill(track.eta());
556  meETLTrackEffPhiMtd_[1]->Fill(track.phi());
557  meETLTrackEffPtMtd_[1]->Fill(track.pt());
558  if (twoETLdiscs) {
559  meETLTrackEffEta2Mtd_[1]->Fill(track.eta());
560  meETLTrackEffPhi2Mtd_[1]->Fill(track.phi());
561  meETLTrackEffPt2Mtd_[1]->Fill(track.pt());
562  }
563  }
564  }
565  }
566 
567  LogDebug("MtdTracksValidation") << "Track p/pt = " << track.p() << " " << track.pt() << " eta " << track.eta()
568  << " BTL " << isBTL << " ETL " << isETL << " 2disks " << twoETLdiscs;
569 
570  // TrackingParticle based matching
571 
572  const reco::TrackBaseRef trkrefb(trackref);
573  auto tp_info = getMatchedTP(trkrefb);
574 
575  meTrackPtTot_->Fill(trackGen.pt());
576  meTrackEtaTot_->Fill(std::abs(trackGen.eta()));
577  if (tp_info != nullptr && mvaTPSel(**tp_info)) {
578  if (track.pt() < 12.) {
579  if (isBTL) {
580  meBTLTrackMatchedTPPtResMtd_->Fill(std::abs(track.pt() - (*tp_info)->pt()) /
581  std::abs(trackGen.pt() - (*tp_info)->pt()));
582  meBTLTrackMatchedTPPtRatioGen_->Fill(trackGen.pt() / (*tp_info)->pt());
583  meBTLTrackMatchedTPPtRatioMtd_->Fill(track.pt() / (*tp_info)->pt());
585  (*tp_info)->pt(), std::abs(track.pt() - (*tp_info)->pt()) / std::abs(trackGen.pt() - (*tp_info)->pt()));
586  meBTLTrackMatchedTPDPtvsPtGen_->Fill((*tp_info)->pt(),
587  (trackGen.pt() - (*tp_info)->pt()) / (*tp_info)->pt());
588  meBTLTrackMatchedTPDPtvsPtMtd_->Fill((*tp_info)->pt(), (track.pt() - (*tp_info)->pt()) / (*tp_info)->pt());
589  }
590  if (isETL && !twoETLdiscs) {
591  meETLTrackMatchedTPPtResMtd_->Fill(std::abs(track.pt() - (*tp_info)->pt()) /
592  std::abs(trackGen.pt() - (*tp_info)->pt()));
593  meETLTrackMatchedTPPtRatioGen_->Fill(trackGen.pt() / (*tp_info)->pt());
594  meETLTrackMatchedTPPtRatioMtd_->Fill(track.pt() / (*tp_info)->pt());
596  (*tp_info)->pt(), std::abs(track.pt() - (*tp_info)->pt()) / std::abs(trackGen.pt() - (*tp_info)->pt()));
597  meETLTrackMatchedTPDPtvsPtGen_->Fill((*tp_info)->pt(),
598  (trackGen.pt() - (*tp_info)->pt()) / ((*tp_info)->pt()));
599  meETLTrackMatchedTPDPtvsPtMtd_->Fill((*tp_info)->pt(),
600  (track.pt() - (*tp_info)->pt()) / ((*tp_info)->pt()));
601  }
602  if (isETL && twoETLdiscs) {
603  meETLTrackMatchedTP2PtResMtd_->Fill(std::abs(track.pt() - (*tp_info)->pt()) /
604  std::abs(trackGen.pt() - (*tp_info)->pt()));
605  meETLTrackMatchedTP2PtRatioGen_->Fill(trackGen.pt() / (*tp_info)->pt());
606  meETLTrackMatchedTP2PtRatioMtd_->Fill(track.pt() / (*tp_info)->pt());
608  (*tp_info)->pt(), std::abs(track.pt() - (*tp_info)->pt()) / std::abs(trackGen.pt() - (*tp_info)->pt()));
609  meETLTrackMatchedTP2DPtvsPtGen_->Fill((*tp_info)->pt(),
610  (trackGen.pt() - (*tp_info)->pt()) / ((*tp_info)->pt()));
611  meETLTrackMatchedTP2DPtvsPtMtd_->Fill((*tp_info)->pt(),
612  (track.pt() - (*tp_info)->pt()) / ((*tp_info)->pt()));
613  }
614  }
615  const bool withMTD = (m_tp2detid.find(*tp_info) != m_tp2detid.end());
616  LogDebug("MtdTracksValidation") << "Matched with selected TP, MTD sim hits association: " << withMTD;
617  if (noCrack) {
618  meTrackMatchedTPEffPtTot_->Fill(trackGen.pt());
619  if (withMTD) {
620  meTrackMatchedTPmtdEffPtTot_->Fill(trackGen.pt());
621  }
622  }
623  meTrackMatchedTPEffEtaTot_->Fill(std::abs(trackGen.eta()));
624  if (withMTD) {
625  meTrackMatchedTPmtdEffEtaTot_->Fill(std::abs(trackGen.eta()));
626  }
627  if (isBTL || isETL) {
628  if (noCrack) {
629  meTrackMatchedTPEffPtMtd_->Fill(trackGen.pt());
630  if (isBTL || twoETLdiscs) {
631  meTrackMatchedTPEffPtEtl2Mtd_->Fill(trackGen.pt());
632  }
633  if (withMTD) {
634  meTrackMatchedTPmtdEffPtMtd_->Fill(trackGen.pt());
635  }
636  }
637  meTrackMatchedTPEffEtaMtd_->Fill(std::abs(trackGen.eta()));
638  if (isBTL || twoETLdiscs) {
639  meTrackMatchedTPEffEtaEtl2Mtd_->Fill(std::abs(trackGen.eta()));
640  }
641  if (withMTD) {
642  meTrackMatchedTPmtdEffEtaMtd_->Fill(std::abs(trackGen.eta()));
643  }
644  }
645 
646  size_t nlayers(0);
647  float extrho(0.);
648  float exteta(0.);
649  float extphi(0.);
650  float selvar(0.);
651  auto accept = checkAcceptance(trackGen, iEvent, iSetup, nlayers, extrho, exteta, extphi, selvar);
652  if (accept.first && std::abs(exteta) < trackMaxBtlEta_) {
654  meExtraBTLeneInCone_->Fill(selvar);
655  }
656  if (accept.second) {
657  if (std::abs(exteta) < trackMaxBtlEta_) {
659  }
660  if (noCrack) {
661  meExtraPtMtd_->Fill(trackGen.pt());
662  if (nlayers == 2) {
663  meExtraPtEtl2Mtd_->Fill(trackGen.pt());
664  }
665  }
666  meExtraEtaMtd_->Fill(std::abs(trackGen.eta()));
667  if (nlayers == 2) {
668  meExtraEtaEtl2Mtd_->Fill(std::abs(trackGen.eta()));
669  }
670  if (accept.first && accept.second && !(isBTL || isETL)) {
671  edm::LogInfo("MtdTracksValidation")
672  << "MtdTracksValidation: extender fail in " << iEvent.id().run() << " " << iEvent.id().event()
673  << " pt= " << trackGen.pt() << " eta= " << trackGen.eta();
674  meExtraMTDfailExtenderEta_->Fill(std::abs(trackGen.eta()));
675  if (noCrack) {
676  meExtraMTDfailExtenderPt_->Fill(trackGen.pt());
677  }
678  }
679  }
680 
681  } // TP matching
682  }
683 
684  if (isGoodVtx) {
685  const bool vtxFake = primRecoVtx.isFake();
686 
687  if (mvaRecSel(trackGen, primRecoVtx, t0Safe[trackref], Sigmat0Safe[trackref])) {
688  // reco-gen matching used for MVA quality flag
689 
690  if (noCrack) {
691  meMVATrackEffPtTot_->Fill(trackGen.pt());
692  }
693  meMVATrackEffEtaTot_->Fill(std::abs(trackGen.eta()));
694 
695  double dZ = trackGen.vz() - zsim;
696  double dT(-9999.);
697  double pullT(-9999.);
698  if (Sigmat0Safe[trackref] != -1.) {
699  dT = t0Safe[trackref] - tsim;
700  pullT = dT / Sigmat0Safe[trackref];
701  }
702  for (const auto& genP : mc->particle_range()) {
703  // select status 1 genParticles and match them to the reconstructed track
704 
705  float charge = pdTable->particle(HepPDT::ParticleID(genP->pdg_id())) != nullptr
706  ? pdTable->particle(HepPDT::ParticleID(genP->pdg_id()))->charge()
707  : 0.f;
708  if (mvaGenSel(*genP, charge)) {
709  if (mvaGenRecMatch(*genP, zsim, trackGen, vtxFake)) {
711  if (noCrack) {
712  meMVATrackMatchedEffPtTot_->Fill(trackGen.pt());
713  }
714  meMVATrackMatchedEffEtaTot_->Fill(std::abs(trackGen.eta()));
715  if (isBTL || isETL) {
716  meMVATrackResTot_->Fill(dT);
717  meMVATrackPullTot_->Fill(pullT);
718  if (noCrack) {
719  meMVATrackMatchedEffPtMtd_->Fill(trackGen.pt());
720  }
721  meMVATrackMatchedEffEtaMtd_->Fill(std::abs(trackGen.eta()));
722  }
723  break;
724  }
725  }
726  }
727  }
728  } // MC truth matich analysis for good PV
729  } //RECO tracks loop
730 }
731 
732 const std::pair<bool, bool> MtdTracksValidation::checkAcceptance(const reco::Track& track,
733  const edm::Event& iEvent,
734  edm::EventSetup const& iSetup,
735  size_t& nlayers,
736  float& extrho,
737  float& exteta,
738  float& extphi,
739  float& selvar) {
740  bool isMatched(false);
741  nlayers = 0;
742  extrho = 0.;
743  exteta = -999.;
744  extphi = -999.;
745  selvar = 0.;
746 
747  auto geometryHandle = iSetup.getTransientHandle(mtdgeoToken_);
748  const MTDGeometry* geom = geometryHandle.product();
749  auto topologyHandle = iSetup.getTransientHandle(mtdtopoToken_);
750  const MTDTopology* topology = topologyHandle.product();
751 
752  auto layerHandle = iSetup.getTransientHandle(mtdlayerToken_);
753  const MTDDetLayerGeometry* layerGeo = layerHandle.product();
754 
755  auto magfieldHandle = iSetup.getTransientHandle(magfieldToken_);
756  const MagneticField* mfield = magfieldHandle.product();
757 
758  auto ttrackBuilder = iSetup.getTransientHandle(builderToken_);
759 
760  auto tTrack = ttrackBuilder->build(track);
761  TrajectoryStateOnSurface tsos = tTrack.outermostMeasurementState();
762  float theMaxChi2 = 500.;
763  float theNSigma = 10.;
764  std::unique_ptr<MeasurementEstimator> theEstimator =
765  std::make_unique<Chi2MeasurementEstimator>(theMaxChi2, theNSigma);
767 
768  auto btlRecHitsHandle = makeValid(iEvent.getHandle(btlRecHitsToken_));
769  auto etlRecHitsHandle = makeValid(iEvent.getHandle(etlRecHitsToken_));
770 
771  edm::LogVerbatim("MtdTracksValidation")
772  << "MtdTracksValidation: extrapolating track, pt= " << track.pt() << " eta= " << track.eta();
773 
774  //try BTL
775  bool inBTL = false;
776  float eneSum(0.);
777  const std::vector<const DetLayer*>& layersBTL = layerGeo->allBTLLayers();
778  for (const DetLayer* ilay : layersBTL) {
779  std::pair<bool, TrajectoryStateOnSurface> comp = ilay->compatible(tsos, prop, *theEstimator);
780  if (!comp.first)
781  continue;
782  if (!inBTL) {
783  inBTL = true;
784  extrho = comp.second.globalPosition().perp();
785  exteta = comp.second.globalPosition().eta();
786  extphi = comp.second.globalPosition().phi();
787  edm::LogVerbatim("MtdTracksValidation") << "MtdTracksValidation: extrapolation at BTL surface, rho= " << extrho
788  << " eta= " << exteta << " phi= " << extphi;
789  }
790  std::vector<DetLayer::DetWithState> compDets = ilay->compatibleDets(tsos, prop, *theEstimator);
791  for (const auto& detWithState : compDets) {
792  const auto& det = detWithState.first;
793 
794  // loop on compatible rechits and check energy in a fixed size cone around the extrapolation point
795 
796  edm::LogVerbatim("MtdTracksValidation")
797  << "MtdTracksValidation: DetId= " << det->geographicalId().rawId()
798  << " gp= " << detWithState.second.globalPosition().x() << " " << detWithState.second.globalPosition().y()
799  << " " << detWithState.second.globalPosition().z() << " rho= " << detWithState.second.globalPosition().perp()
800  << " eta= " << detWithState.second.globalPosition().eta()
801  << " phi= " << detWithState.second.globalPosition().phi();
802 
803  for (const auto& recHit : *btlRecHitsHandle) {
804  BTLDetId detId = recHit.id();
805  DetId geoId = detId.geographicalId(MTDTopologyMode::crysLayoutFromTopoMode(topology->getMTDTopologyMode()));
806  const MTDGeomDet* thedet = geom->idToDet(geoId);
807  if (thedet == nullptr)
808  throw cms::Exception("MtdTracksValidation") << "GeographicalID: " << std::hex << geoId.rawId() << " ("
809  << detId.rawId() << ") is invalid!" << std::dec << std::endl;
810  if (geoId == det->geographicalId()) {
811  const ProxyMTDTopology& topoproxy = static_cast<const ProxyMTDTopology&>(thedet->topology());
812  const RectangularMTDTopology& topo = static_cast<const RectangularMTDTopology&>(topoproxy.specificTopology());
813 
814  Local3DPoint local_point(0., 0., 0.);
815  local_point = topo.pixelToModuleLocalPoint(local_point, detId.row(topo.nrows()), detId.column(topo.nrows()));
816  const auto& global_point = thedet->toGlobal(local_point);
817  edm::LogVerbatim("MtdTracksValidation")
818  << "MtdTracksValidation: Hit id= " << detId.rawId() << " ene= " << recHit.energy()
819  << " dr= " << reco::deltaR(global_point, detWithState.second.globalPosition());
820  if (reco::deltaR(global_point, detWithState.second.globalPosition()) < cluDRradius_) {
821  eneSum += recHit.energy();
822  //extrho = detWithState.second.globalPosition().perp();
823  //exteta = detWithState.second.globalPosition().eta();
824  //extphi = detWithState.second.globalPosition().phi();
825  }
826  }
827  }
828  }
829  if (eneSum > depositBTLthreshold_) {
830  nlayers++;
831  selvar = eneSum;
832  isMatched = true;
833  edm::LogVerbatim("MtdTracksValidation")
834  << "MtdTracksValidation: BTL matched, energy= " << eneSum << " #layers= " << nlayers;
835  }
836  }
837  if (inBTL) {
838  return std::make_pair(inBTL, isMatched);
839  }
840 
841  //try ETL
842  bool inETL = false;
843  const std::vector<const DetLayer*>& layersETL = layerGeo->allETLLayers();
844  for (const DetLayer* ilay : layersETL) {
845  size_t hcount(0);
846  const BoundDisk& disk = static_cast<const MTDSectorForwardDoubleLayer*>(ilay)->specificSurface();
847  const double diskZ = disk.position().z();
848  if (tsos.globalPosition().z() * diskZ < 0)
849  continue; // only propagate to the disk that's on the same side
850  std::pair<bool, TrajectoryStateOnSurface> comp = ilay->compatible(tsos, prop, *theEstimator);
851  if (!comp.first)
852  continue;
853  if (!inETL) {
854  inETL = true;
855  extrho = comp.second.globalPosition().perp();
856  exteta = comp.second.globalPosition().eta();
857  extphi = comp.second.globalPosition().phi();
858  }
859  edm::LogVerbatim("MtdTracksValidation") << "MtdTracksValidation: extrapolation at ETL surface, rho= " << extrho
860  << " eta= " << exteta << " phi= " << extphi;
861  std::vector<DetLayer::DetWithState> compDets = ilay->compatibleDets(tsos, prop, *theEstimator);
862  for (const auto& detWithState : compDets) {
863  const auto& det = detWithState.first;
864 
865  // loop on compatible rechits and check hits in a fixed size cone around the extrapolation point
866 
867  edm::LogVerbatim("MtdTracksValidation")
868  << "MtdTracksValidation: DetId= " << det->geographicalId().rawId()
869  << " gp= " << detWithState.second.globalPosition().x() << " " << detWithState.second.globalPosition().y()
870  << " " << detWithState.second.globalPosition().z() << " rho= " << detWithState.second.globalPosition().perp()
871  << " eta= " << detWithState.second.globalPosition().eta()
872  << " phi= " << detWithState.second.globalPosition().phi();
873 
874  for (const auto& recHit : *etlRecHitsHandle) {
875  ETLDetId detId = recHit.id();
876  DetId geoId = detId.geographicalId();
877  const MTDGeomDet* thedet = geom->idToDet(geoId);
878  if (thedet == nullptr)
879  throw cms::Exception("MtdTracksValidation") << "GeographicalID: " << std::hex << geoId.rawId() << " ("
880  << detId.rawId() << ") is invalid!" << std::dec << std::endl;
881  if (geoId == det->geographicalId()) {
882  const ProxyMTDTopology& topoproxy = static_cast<const ProxyMTDTopology&>(thedet->topology());
883  const RectangularMTDTopology& topo = static_cast<const RectangularMTDTopology&>(topoproxy.specificTopology());
884 
885  Local3DPoint local_point(topo.localX(recHit.row()), topo.localY(recHit.column()), 0.);
886  const auto& global_point = thedet->toGlobal(local_point);
887  edm::LogVerbatim("MtdTracksValidation")
888  << "MtdTracksValidation: Hit id= " << detId.rawId() << " time= " << recHit.time()
889  << " dr= " << reco::deltaR(global_point, detWithState.second.globalPosition());
890  if (reco::deltaR(global_point, detWithState.second.globalPosition()) < cluDRradius_) {
891  hcount++;
892  if (hcount == 1) {
893  //extrho = detWithState.second.globalPosition().perp();
894  //exteta = detWithState.second.globalPosition().eta();
895  //extphi = detWithState.second.globalPosition().phi();
896  }
897  }
898  }
899  }
900  }
901  if (hcount > 0) {
902  nlayers++;
903  selvar = (float)hcount;
904  isMatched = true;
905  edm::LogVerbatim("MtdTracksValidation")
906  << "MtdTracksValidation: ETL matched, counts= " << hcount << " #layers= " << nlayers;
907  }
908  }
909 
910  if (!inBTL && !inETL) {
911  edm::LogVerbatim("MtdTracksValidation")
912  << "MtdTracksValidation: track not extrapolating to MTD: pt= " << track.pt() << " eta= " << track.eta()
913  << " phi= " << track.phi() << " vz= " << track.vz()
914  << " vxy= " << std::sqrt(track.vx() * track.vx() + track.vy() * track.vy());
915  }
916  return std::make_pair(inETL, isMatched);
917 }
918 
919 // ------------ method for histogram booking ------------
921  ibook.setCurrentFolder(folder_);
922 
923  // histogram booking
924  meBTLTrackRPTime_ = ibook.book1D("TrackBTLRPTime", "Track t0 with respect to R.P.;t0 [ns]", 100, -1, 3);
925  meBTLTrackEffEtaTot_ = ibook.book1D("TrackBTLEffEtaTot", "Track efficiency vs eta (Tot);#eta_{RECO}", 100, -1.6, 1.6);
927  ibook.book1D("TrackBTLEffPhiTot", "Track efficiency vs phi (Tot);#phi_{RECO} [rad]", 100, -3.2, 3.2);
928  meBTLTrackEffPtTot_ = ibook.book1D("TrackBTLEffPtTot", "Track efficiency vs pt (Tot);pt_{RECO} [GeV]", 50, 0, 10);
929  meBTLTrackEffEtaMtd_ = ibook.book1D("TrackBTLEffEtaMtd", "Track efficiency vs eta (Mtd);#eta_{RECO}", 100, -1.6, 1.6);
931  ibook.book1D("TrackBTLEffPhiMtd", "Track efficiency vs phi (Mtd);#phi_{RECO} [rad]", 100, -3.2, 3.2);
932  meBTLTrackEffPtMtd_ = ibook.book1D("TrackBTLEffPtMtd", "Track efficiency vs pt (Mtd);pt_{RECO} [GeV]", 50, 0, 10);
934  ibook.book1D("TrackBTLPtRes", "Track pT resolution ;pT_{Gentrack}-pT_{MTDtrack}/pT_{Gentrack} ", 100, -0.1, 0.1);
935  meETLTrackRPTime_ = ibook.book1D("TrackETLRPTime", "Track t0 with respect to R.P.;t0 [ns]", 100, -1, 3);
937  ibook.book1D("TrackETLEffEtaTotZneg", "Track efficiency vs eta (Tot) (-Z);#eta_{RECO}", 100, -3.2, -1.4);
939  ibook.book1D("TrackETLEffEtaTotZpos", "Track efficiency vs eta (Tot) (+Z);#eta_{RECO}", 100, 1.4, 3.2);
941  ibook.book1D("TrackETLEffPhiTotZneg", "Track efficiency vs phi (Tot) (-Z);#phi_{RECO} [rad]", 100, -3.2, 3.2);
943  ibook.book1D("TrackETLEffPhiTotZpos", "Track efficiency vs phi (Tot) (+Z);#phi_{RECO} [rad]", 100, -3.2, 3.2);
945  ibook.book1D("TrackETLEffPtTotZneg", "Track efficiency vs pt (Tot) (-Z);pt_{RECO} [GeV]", 50, 0, 10);
947  ibook.book1D("TrackETLEffPtTotZpos", "Track efficiency vs pt (Tot) (+Z);pt_{RECO} [GeV]", 50, 0, 10);
949  ibook.book1D("TrackETLEffEtaMtdZneg", "Track efficiency vs eta (Mtd) (-Z);#eta_{RECO}", 100, -3.2, -1.4);
951  ibook.book1D("TrackETLEffEtaMtdZpos", "Track efficiency vs eta (Mtd) (+Z);#eta_{RECO}", 100, 1.4, 3.2);
953  ibook.book1D("TrackETLEffPhiMtdZneg", "Track efficiency vs phi (Mtd) (-Z);#phi_{RECO} [rad]", 100, -3.2, 3.2);
955  ibook.book1D("TrackETLEffPhiMtdZpos", "Track efficiency vs phi (Mtd) (+Z);#phi_{RECO} [rad]", 100, -3.2, 3.2);
957  ibook.book1D("TrackETLEffPtMtdZneg", "Track efficiency vs pt (Mtd) (-Z);pt_{RECO} [GeV]", 50, 0, 10);
959  ibook.book1D("TrackETLEffPtMtdZpos", "Track efficiency vs pt (Mtd) (+Z);pt_{RECO} [GeV]", 50, 0, 10);
961  ibook.book1D("TrackETLEffEta2MtdZneg", "Track efficiency vs eta (Mtd 2 hit) (-Z);#eta_{RECO}", 100, -3.2, -1.4);
963  ibook.book1D("TrackETLEffEta2MtdZpos", "Track efficiency vs eta (Mtd 2 hit) (+Z);#eta_{RECO}", 100, 1.4, 3.2);
964  meETLTrackEffPhi2Mtd_[0] = ibook.book1D(
965  "TrackETLEffPhi2MtdZneg", "Track efficiency vs phi (Mtd 2 hit) (-Z);#phi_{RECO} [rad]", 100, -3.2, 3.2);
966  meETLTrackEffPhi2Mtd_[1] = ibook.book1D(
967  "TrackETLEffPhi2MtdZpos", "Track efficiency vs phi (Mtd 2 hit) (+Z);#phi_{RECO} [rad]", 100, -3.2, 3.2);
969  ibook.book1D("TrackETLEffPt2MtdZneg", "Track efficiency vs pt (Mtd 2 hit) (-Z);pt_{RECO} [GeV]", 50, 0, 10);
971  ibook.book1D("TrackETLEffPt2MtdZpos", "Track efficiency vs pt (Mtd 2 hit) (+Z);pt_{RECO} [GeV]", 50, 0, 10);
973  ibook.book1D("TrackETLPtRes", "Track pT resolution;pT_{Gentrack}-pT_{MTDtrack}/pT_{Gentrack} ", 100, -0.1, 0.1);
974 
975  meTracktmtd_ = ibook.book1D("Tracktmtd", "Track time from TrackExtenderWithMTD;tmtd [ns]", 150, 1, 16);
976  meTrackt0Src_ = ibook.book1D("Trackt0Src", "Track time from TrackExtenderWithMTD;t0Src [ns]", 100, -1.5, 1.5);
978  ibook.book1D("TrackSigmat0Src", "Time Error from TrackExtenderWithMTD; #sigma_{t0Src} [ns]", 100, 0, 0.1);
979 
980  meTrackt0Pid_ = ibook.book1D("Trackt0Pid", "Track t0 as stored in TofPid;t0 [ns]", 100, -1, 1);
981  meTrackSigmat0Pid_ = ibook.book1D("TrackSigmat0Pid", "Sigmat0 as stored in TofPid; #sigma_{t0} [ns]", 100, 0, 0.1);
982  meTrackt0SafePid_ = ibook.book1D("Trackt0SafePID", "Track t0 Safe as stored in TofPid;t0 [ns]", 100, -1, 1);
984  ibook.book1D("TrackSigmat0SafePID", "Sigmat0 Safe as stored in TofPid; #sigma_{t0} [ns]", 100, 0, 0.1);
985  meTrackNumHits_ = ibook.book1D("TrackNumHits", "Number of valid MTD hits per track ; Number of hits", 10, -5, 5);
986  meTrackNumHitsNT_ = ibook.book1D(
987  "TrackNumHitsNT", "Number of valid MTD hits per track no time associated; Number of hits", 10, -5, 5);
988  meTrackMVAQual_ = ibook.book1D("TrackMVAQual", "Track MVA Quality as stored in Value Map ; MVAQual", 100, 0, 1);
990  "TrackPathLenghtvsEta", "MTD Track pathlength vs MTD track Eta;|#eta|;Pathlength", 100, 0, 3.2, 100.0, 400.0, "S");
991 
992  meMVATrackEffPtTot_ = ibook.book1D("MVAEffPtTot", "Pt of tracks associated to LV; track pt [GeV] ", 110, 0., 11.);
994  ibook.book1D("MVAMatchedEffPtTot", "Pt of tracks associated to LV matched to GEN; track pt [GeV] ", 110, 0., 11.);
996  "MVAMatchedEffPtMtd", "Pt of tracks associated to LV matched to GEN with time; track pt [GeV] ", 110, 0., 11.);
997 
998  meExtraPtMtd_ = ibook.book1D("ExtraPtMtd", "Pt of tracks extrapolated to hits; track pt [GeV] ", 110, 0., 11.);
1000  ibook.book1D("ExtraPtEtl2Mtd", "Pt of tracks extrapolated to hits, 2 ETL layers; track pt [GeV] ", 110, 0., 11.);
1001 
1002  meTrackPtTot_ = ibook.book1D("TrackPtTot", "Pt of tracks ; track pt [GeV] ", 110, 0., 11.);
1004  ibook.book1D("MatchedTPEffPtTot", "Pt of tracks matched to TP; track pt [GeV] ", 110, 0., 11.);
1006  ibook.book1D("MatchedTPEffPtMtd", "Pt of tracks matched to TP with time; track pt [GeV] ", 110, 0., 11.);
1008  "MatchedTPEffPtEtl2Mtd", "Pt of tracks matched to TP with time, 2 ETL hits; track pt [GeV] ", 110, 0., 11.);
1009 
1011  "TrackMatchedTPBTLPtResMtd",
1012  "Pt resolution of tracks matched to TP-BTL hit ;|pT_{MTDtrack}-pT_{truth}|/|pT_{Gentrack}-pT_{truth}| ",
1013  100,
1014  0.,
1015  4.);
1017  "TrackMatchedTPETLPtResMtd",
1018  "Pt resolution of tracks matched to TP-ETL hit ;|pT_{MTDtrack}-pT_{truth}|/|pT_{Gentrack}-pT_{truth}| ",
1019  100,
1020  0.,
1021  4.);
1023  "TrackMatchedTPETL2PtResMtd",
1024  "Pt resolution of tracks matched to TP-ETL 2hits ;|pT_{MTDtrack}-pT_{truth}|/|pT_{Gentrack}-pT_{truth}| ",
1025  100,
1026  0.,
1027  4.);
1029  "TrackMatchedTPBTLPtRatioGen", "Pt ratio of Gentracks (BTL) ;pT_{Gentrack}/pT_{truth} ", 100, 0.9, 1.1);
1031  "TrackMatchedTPETLPtRatioGen", "Pt ratio of Gentracks (ETL 1hit) ;pT_{Gentrack}/pT_{truth} ", 100, 0.9, 1.1);
1033  "TrackMatchedTPETL2PtRatioGen", "Pt ratio of Gentracks (ETL 2hits) ;pT_{Gentrack}/pT_{truth} ", 100, 0.9, 1.1);
1034  meBTLTrackMatchedTPPtRatioMtd_ = ibook.book1D("TrackMatchedTPBTLPtRatioMtd",
1035  "Pt ratio of tracks matched to TP-BTL hits ;pT_{MTDtrack}/pT_{truth} ",
1036  100,
1037  0.9,
1038  1.1);
1039  meETLTrackMatchedTPPtRatioMtd_ = ibook.book1D("TrackMatchedTPETLPtRatioMtd",
1040  "Pt ratio of tracks matched to TP-ETL hits ;pT_{MTDtrack}/pT_{truth} ",
1041  100,
1042  0.9,
1043  1.1);
1045  ibook.book1D("TrackMatchedTPETL2PtRatioMtd",
1046  "Pt ratio of tracks matched to TP-ETL 2hits ;pT_{MTDtrack}/pT_{truth} ",
1047  100,
1048  0.9,
1049  1.1);
1050  meBTLTrackMatchedTPPtResvsPtMtd_ = ibook.bookProfile("TrackMatchedTPBTLPtResvsPtMtd",
1051  "Pt resolution of tracks matched to TP-BTL hit vs Pt;pT_{truth} "
1052  "[GeV];|pT_{MTDtrack}-pT_{truth}|/|pT_{Gentrack}-pT_{truth}| ",
1053  20,
1054  0.7,
1055  10.,
1056  0.,
1057  4.,
1058  "s");
1059  meETLTrackMatchedTPPtResvsPtMtd_ = ibook.bookProfile("TrackMatchedTPETLPtResvsPtMtd",
1060  "Pt resolution of tracks matched to TP-ETL hit vs Pt;pT_{truth} "
1061  "[GeV];|pT_{MTDtrack}-pT_{truth}|/|pT_{Gentrack}-pT_{truth}| ",
1062  20,
1063  0.7,
1064  10.,
1065  0.,
1066  4.,
1067  "s");
1069  ibook.bookProfile("TrackMatchedTPETL2PtResvsPtMtd",
1070  "Pt resolution of tracks matched to TP-ETL 2hits Pt pT;pT_{truth} "
1071  "[GeV];|pT_{MTDtrack}-pT_{truth}|/|pT_{Gentrack}-pT_{truth}| ",
1072  20,
1073  0.7,
1074  10.,
1075  0.,
1076  4.,
1077  "s");
1079  "TrackMatchedTPBTLDPtvsPtGen",
1080  "Pt relative difference of Gentracks (BTL) vs Pt;pT_{truth} [GeV];pT_{Gentrack}-pT_{truth}/pT_{truth} ",
1081  20,
1082  0.7,
1083  10.,
1084  -0.1,
1085  0.1,
1086  "s");
1088  "TrackMatchedTPETLDPtvsPtGen",
1089  "Pt relative difference of Gentracks (ETL 1hit) vs Pt;pT_{truth} [GeV];pT_{Gentrack}-pT_{truth}/pT_{truth} ",
1090  20,
1091  0.7,
1092  10.,
1093  -0.1,
1094  0.1,
1095  "s");
1097  "TrackMatchedTPETL2DPtvsPtGen",
1098  "Pt relative difference of Gentracks (ETL 2hits) vs Pt;pT_{truth} [GeV];pT_{Gentrack}-pT_{truth}/pT_{truth} ",
1099  20,
1100  0.7,
1101  10.,
1102  -0.1,
1103  0.1,
1104  "s");
1105  meBTLTrackMatchedTPDPtvsPtMtd_ = ibook.bookProfile("TrackMatchedTPBTLDPtvsPtMtd",
1106  "Pt relative difference of tracks matched to TP-BTL hits vs "
1107  "Pt;pT_{truth} [GeV];pT_{MTDtrack}-pT_{truth}/pT_{truth} ",
1108  20,
1109  0.7,
1110  10.,
1111  -0.1,
1112  0.1,
1113  "s");
1114  meETLTrackMatchedTPDPtvsPtMtd_ = ibook.bookProfile("TrackMatchedTPETLDPtvsPtMtd",
1115  "Pt relative difference of tracks matched to TP-ETL hits vs "
1116  "Pt;pT_{truth} [GeV];pT_{MTDtrack}-pT_{truth}/pT_{truth} ",
1117  20,
1118  0.7,
1119  10.,
1120  -0.1,
1121  0.1,
1122  "s");
1123  meETLTrackMatchedTP2DPtvsPtMtd_ = ibook.bookProfile("TrackMatchedTPETL2DPtvsPtMtd",
1124  "Pt relative difference of tracks matched to TP-ETL 2hits vs "
1125  "Pt;pT_{truth} [GeV];pT_{MTDtrack}-pT_{truth}/pT_{truth} ",
1126  20,
1127  0.7,
1128  10.,
1129  -0.1,
1130  0.1,
1131  "s");
1132 
1134  ibook.book1D("MatchedTPmtdEffPtTot", "Pt of tracks matched to TP-mtd hit; track pt [GeV] ", 110, 0., 11.);
1136  "MatchedTPmtdEffPtMtd", "Pt of tracks matched to TP-mtd hit with time; track pt [GeV] ", 110, 0., 11.);
1137 
1138  meMVATrackEffEtaTot_ = ibook.book1D("MVAEffEtaTot", "Eta of tracks associated to LV; track eta ", 66, 0., 3.3);
1140  ibook.book1D("MVAMatchedEffEtaTot", "Eta of tracks associated to LV matched to GEN; track eta ", 66, 0., 3.3);
1142  "MVAMatchedEffEtaMtd", "Eta of tracks associated to LV matched to GEN with time; track eta ", 66, 0., 3.3);
1143 
1144  meExtraEtaMtd_ = ibook.book1D("ExtraEtaMtd", "Eta of tracks extrapolated to hits; track eta ", 66, 0., 3.3);
1146  ibook.book1D("ExtraEtaEtl2Mtd", "Eta of tracks extrapolated to hits, 2 ETL layers; track eta ", 66, 0., 3.3);
1147 
1148  meTrackEtaTot_ = ibook.book1D("TrackEtaTot", "Eta of tracks ; track eta ", 66, 0., 3.3);
1150  ibook.book1D("MatchedTPEffEtaTot", "Eta of tracks matched to TP; track eta ", 66, 0., 3.3);
1151  meMVATrackEffEtaTot_ = ibook.book1D("MVAEffEtaTot", "Eta of tracks ; track eta ", 66, 0., 3.3);
1153  ibook.book1D("MatchedTPEffEtaMtd", "Eta of tracks matched to TP with time; track eta ", 66, 0., 3.3);
1155  "MatchedTPEffEtaEtl2Mtd", "Eta of tracks matched to TP with time, 2 ETL hits; track eta ", 66, 0., 3.3);
1156 
1158  ibook.book1D("MatchedTPmtdEffEtaTot", "Eta of tracks matched to TP-mtd hit; track eta ", 66, 0., 3.3);
1160  ibook.book1D("MatchedTPmtdEffEtaMtd", "Eta of tracks matched to TP-mtd hit with time; track eta ", 66, 0., 3.3);
1161 
1162  meMVATrackResTot_ = ibook.book1D(
1163  "MVATrackRes", "t_{rec} - t_{sim} for LV associated tracks; t_{rec} - t_{sim} [ns] ", 120, -0.15, 0.15);
1165  ibook.book1D("MVATrackPull", "Pull for associated tracks; (t_{rec}-t_{sim})/#sigma_{t}", 50, -5., 5.);
1166  meMVATrackZposResTot_ = ibook.book1D(
1167  "MVATrackZposResTot", "Z_{PCA} - Z_{sim} for associated tracks;Z_{PCA} - Z_{sim} [cm] ", 100, -0.1, 0.1);
1168 
1170  ibook.book1D("ExtraPhiAtBTL", "Phi at BTL surface of extrapolated tracks; phi [deg]", 720, -180., 180.);
1171  meExtraPhiAtBTLmatched_ = ibook.book1D("ExtraPhiAtBTLmatched",
1172  "Phi at BTL surface of extrapolated tracksi matched with BTL hits; phi [deg]",
1173  720,
1174  -180.,
1175  180.);
1176  meExtraBTLeneInCone_ = ibook.book1D(
1177  "ExtraBTLeneInCone", "BTL reconstructed energy in cone arounnd extrapolated track; E [MeV]", 100, 0., 50.);
1179  ibook.book1D("ExtraMTDfailExtenderEta",
1180  "Eta of tracks extrapolated to MTD with no track extender match to hits; track eta",
1181  66,
1182  0.,
1183  3.3);
1184  ;
1186  ibook.book1D("ExtraMTDfailExtenderPt",
1187  "Pt of tracks extrapolated to MTD with no track extender match to hits; track pt [GeV] ",
1188  110,
1189  0.,
1190  11.);
1191 }
1192 
1193 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
1194 
1197 
1198  desc.add<std::string>("folder", "MTD/Tracks");
1199  desc.add<edm::InputTag>("inputTagG", edm::InputTag("generalTracks"));
1200  desc.add<edm::InputTag>("inputTagT", edm::InputTag("trackExtenderWithMTD"));
1201  desc.add<edm::InputTag>("inputTagV", edm::InputTag("offlinePrimaryVertices4D"));
1202  desc.add<edm::InputTag>("inputTagH", edm::InputTag("generatorSmeared"));
1203  desc.add<edm::InputTag>("SimTag", edm::InputTag("mix", "MergedTrackTruth"));
1204  desc.add<edm::InputTag>("TPtoRecoTrackAssoc", edm::InputTag("trackingParticleRecoTrackAsssociation"));
1205  desc.add<edm::InputTag>("btlSimHits", edm::InputTag("mix", "g4SimHitsFastTimerHitsBarrel"));
1206  desc.add<edm::InputTag>("etlSimHits", edm::InputTag("mix", "g4SimHitsFastTimerHitsEndcap"));
1207  desc.add<edm::InputTag>("btlRecHits", edm::InputTag("mtdRecHits", "FTLBarrel"));
1208  desc.add<edm::InputTag>("etlRecHits", edm::InputTag("mtdRecHits", "FTLEndcap"));
1209  desc.add<edm::InputTag>("tmtd", edm::InputTag("trackExtenderWithMTD:generalTracktmtd"));
1210  desc.add<edm::InputTag>("sigmatmtd", edm::InputTag("trackExtenderWithMTD:generalTracksigmatmtd"));
1211  desc.add<edm::InputTag>("t0Src", edm::InputTag("trackExtenderWithMTD:generalTrackt0"));
1212  desc.add<edm::InputTag>("sigmat0Src", edm::InputTag("trackExtenderWithMTD:generalTracksigmat0"));
1213  desc.add<edm::InputTag>("trackAssocSrc", edm::InputTag("trackExtenderWithMTD:generalTrackassoc"))
1214  ->setComment("Association between General and MTD Extended tracks");
1215  desc.add<edm::InputTag>("pathLengthSrc", edm::InputTag("trackExtenderWithMTD:generalTrackPathLength"));
1216  desc.add<edm::InputTag>("t0SafePID", edm::InputTag("tofPID:t0safe"));
1217  desc.add<edm::InputTag>("sigmat0SafePID", edm::InputTag("tofPID:sigmat0safe"));
1218  desc.add<edm::InputTag>("sigmat0PID", edm::InputTag("tofPID:sigmat0"));
1219  desc.add<edm::InputTag>("t0PID", edm::InputTag("tofPID:t0"));
1220  desc.add<edm::InputTag>("trackMVAQual", edm::InputTag("mtdTrackQualityMVA:mtdQualMVA"));
1221  desc.add<double>("trackMinimumPt", 0.7); // [GeV]
1222  desc.add<double>("trackMaximumBtlEta", 1.5);
1223  desc.add<double>("trackMinimumEtlEta", 1.6);
1224  desc.add<double>("trackMaximumEtlEta", 3.);
1225  desc.addUntracked<bool>("optionalPlots", true);
1226 
1227  descriptions.add("mtdTracksValid", desc);
1228 }
1229 
1231  bool match = false;
1232  if (gp.status() != 1) {
1233  return match;
1234  }
1235  match = charge != 0.f && gp.momentum().perp() > pTcut_ && std::abs(gp.momentum().eta()) < etacutGEN_;
1236  return match;
1237 }
1238 
1240  bool match = false;
1241  if (tp.status() != 1) {
1242  return match;
1243  }
1244  auto x_pv = tp.parentVertex()->position().x();
1245  auto y_pv = tp.parentVertex()->position().y();
1246  auto z_pv = tp.parentVertex()->position().z();
1247 
1248  auto r_pv = std::sqrt(x_pv * x_pv + y_pv * y_pv);
1249 
1250  match = tp.charge() != 0 && tp.pt() > pTcut_ && std::abs(tp.eta()) < etacutGEN_ && r_pv < rBTL_ && z_pv < zETL_;
1251  return match;
1252 }
1253 
1255  const reco::Vertex& vtx,
1256  const double& t0,
1257  const double& st0) {
1258  bool match = false;
1259  match = trk.pt() > pTcut_ && std::abs(trk.eta()) < etacutREC_ &&
1260  (std::abs(trk.vz() - vtx.z()) <= deltaZcut_ || vtx.isFake());
1261  if (st0 > 0.) {
1262  match = match && std::abs(t0 - vtx.t()) < 3. * st0;
1263  }
1264  return match;
1265 }
1266 
1268  const double& zsim,
1269  const reco::TrackBase& trk,
1270  const bool& vtxFake) {
1271  bool match = false;
1272  double dR = reco::deltaR(genP.momentum(), trk.momentum());
1273  double genPT = genP.momentum().perp();
1274  match = std::abs(genPT - trk.pt()) < trk.pt() * deltaPTcut_ && dR < deltaDRcut_ &&
1275  (std::abs(trk.vz() - zsim) < deltaZcut_ || vtxFake);
1276  return match;
1277 }
1278 
1280  auto found = r2s_->find(recoTrack);
1281 
1282  // reco track not matched to any TP
1283  if (found == r2s_->end())
1284  return nullptr;
1285 
1286  //matched TP equal to any TP associated to in time events
1287  for (const auto& tp : found->val) {
1288  if (tp.first->eventId().bunchCrossing() == 0)
1289  return &tp.first;
1290  }
1291 
1292  // reco track not matched to any TP from vertex
1293  return nullptr;
1294 }
1295 
MonitorElement * meETLTrackEffEtaTot_[2]
edm::EDGetTokenT< edm::ValueMap< int > > trackAssocToken_
static constexpr double etacutGEN_
Log< level::Info, true > LogVerbatim
MonitorElement * meBTLTrackMatchedTPPtResMtd_
MonitorElement * meTrackMatchedTPEffEtaMtd_
MonitorElement * meBTLTrackEffPhiMtd_
edm::EDGetTokenT< edm::ValueMap< float > > trackMVAQualToken_
edm::EDGetTokenT< edm::ValueMap< float > > pathLengthToken_
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
MonitorElement * meMVATrackResTot_
MonitorElement * meETLTrackMatchedTPPtResvsPtMtd_
MonitorElement * meTrackNumHits_
MonitorElement * meTrackt0Pid_
MonitorElement * meExtraMTDfailExtenderEta_
static constexpr double pTcut_
const std::pair< bool, bool > checkAcceptance(const reco::Track &, const edm::Event &, const edm::EventSetup &, size_t &, float &, float &, float &, float &)
MonitorElement * meETLTrackMatchedTP2DPtvsPtGen_
MonitorElement * meETLTrackRPTime_
MonitorElement * meETLTrackMatchedTP2PtResMtd_
edm::EDGetTokenT< edm::ValueMap< float > > t0SafePidToken_
std::string folder_
edm::EDGetTokenT< edm::ValueMap< float > > t0SrcToken_
static constexpr double rBTL_
MonitorElement * meETLTrackMatchedTP2PtRatioMtd_
edm::ESGetToken< MTDTopology, MTDTopologyRcd > mtdtopoToken_
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
const std::string folder_
edm::EDGetTokenT< edm::HepMCProduct > HepMCProductToken_
HepPDT::ParticleDataTable ParticleDataTable
T z() const
Definition: PV3DBase.h:61
MonitorElement * meMVATrackMatchedEffEtaMtd_
const bool mvaRecSel(const reco::TrackBase &, const reco::Vertex &, const double &, const double &)
virtual const Topology & topology() const
Definition: GeomDet.cc:67
edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > builderToken_
constexpr NumType convertRadToDeg(NumType radians)
Definition: angle_units.h:21
virtual const PixelTopology & specificTopology() const
MonitorElement * meETLTrackMatchedTP2PtResvsPtMtd_
static constexpr double etacutREC_
static constexpr double deltaPTcut_
MonitorElement * meBTLTrackPtRes_
edm::EDGetTokenT< std::vector< reco::Vertex > > RecVertexToken_
MonitorElement * meETLTrackEffPtTot_[2]
MonitorElement * meBTLTrackEffEtaMtd_
MonitorElement * meETLTrackEffPhi2Mtd_[2]
Definition: MTDHit.h:4
edm::EDGetTokenT< edm::ValueMap< float > > tmtdToken_
MonitorElement * meExtraBTLeneInCone_
MonitorElement * meExtraMTDfailExtenderPt_
MonitorElement * meMVATrackMatchedEffPtTot_
MonitorElement * meMVATrackEffPtTot_
MonitorElement * meExtraEtaMtd_
MonitorElement * meMVATrackEffEtaTot_
const reco::SimToRecoCollection * s2r_
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:31
const unsigned long int uniqueId(const uint32_t x, const EncodedEventId &y)
MonitorElement * meBTLTrackMatchedTPPtRatioMtd_
LocalPoint pixelToModuleLocalPoint(const LocalPoint &plp, int row, int col) const
MonitorElement * meTrackMatchedTPmtdEffEtaTot_
static constexpr double cluDRradius_
const reco::RecoToSimCollection * r2s_
const edm::Ref< std::vector< TrackingParticle > > * getMatchedTP(const reco::TrackBaseRef &)
const_iterator find(const key_type &k) const
find element with specified reference key
const bool mvaGenSel(const HepMC::GenParticle &, const float &)
Detector identifier base class for the MIP Timing Layer.
Definition: MTDDetId.h:21
const bool mvaTPSel(const TrackingParticle &)
MonitorElement * meETLTrackMatchedTP2PtRatioGen_
constexpr NumType convertUnitsTo(double desiredUnits, NumType val)
Definition: GeantUnits.h:73
edm::EDGetTokenT< reco::TrackCollection > GenRecTrackToken_
MonitorElement * meETLTrackEffPhiMtd_[2]
const_iterator end() const
last iterator over the map (read only)
void Fill(long long x)
double pt() const
track transverse momentum
Definition: TrackBase.h:637
edm::EDGetTokenT< reco::TrackCollection > RecTrackToken_
MonitorElement * meExtraEtaEtl2Mtd_
MonitorElement * meBTLTrackEffPhiTot_
edm::ESGetToken< MTDDetLayerGeometry, MTDRecoGeometryRecord > mtdlayerToken_
int iEvent
Definition: GenABIO.cc:224
GlobalPoint globalPosition() const
MonitorElement * meTrackMatchedTPmtdEffEtaMtd_
MonitorElement * meTrackMatchedTPEffEtaTot_
float localX(const float mpX) const override
edm::EDGetTokenT< FTLRecHitCollection > etlRecHitsToken_
MonitorElement * meBTLTrackRPTime_
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:661
MonitorElement * meETLTrackEffEta2Mtd_[2]
MonitorElement * meBTLTrackMatchedTPDPtvsPtMtd_
static constexpr double depositBTLthreshold_
MonitorElement * meETLTrackPtRes_
MonitorElement * meETLTrackEffPt2Mtd_[2]
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:408
MonitorElement * meETLTrackMatchedTPDPtvsPtGen_
MonitorElement * meTrackt0SafePid_
T sqrt(T t)
Definition: SSEVec.h:19
MonitorElement * meBTLTrackEffEtaTot_
MonitorElement * meETLTrackMatchedTPDPtvsPtMtd_
static constexpr double deltaZcut_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isETL(const double eta) const
bool isMatched(TrackingRecHit const &hit)
MonitorElement * meBTLTrackMatchedTPPtResvsPtMtd_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
int nrows() const override
MonitorElement * meMVATrackMatchedEffEtaTot_
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
edm::EDGetTokenT< edm::ValueMap< float > > t0PidToken_
MonitorElement * meTrackMatchedTPEffPtEtl2Mtd_
MonitorElement * meTrackMatchedTPmtdEffPtMtd_
edm::EDGetTokenT< FTLRecHitCollection > btlRecHitsToken_
MonitorElement * meETLTrackMatchedTPPtRatioMtd_
MonitorElement * meTrackMVAQual_
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
MonitorElement * meETLTrackEffPhiTot_[2]
float localY(const float mpY) const override
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
const std::vector< const DetLayer * > & allBTLLayers() const
return the BTL DetLayers (barrel), inside-out
edm::EDGetTokenT< edm::ValueMap< float > > SigmatmtdToken_
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
edm::EDGetTokenT< TrackingParticleCollection > trackingParticleCollectionToken_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
MonitorElement * meTrackNumHitsNT_
Log< level::Info, false > LogInfo
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:652
Definition: DetId.h:17
constexpr NumType convertMmToCm(NumType millimeters)
Definition: angle_units.h:44
edm::EDGetTokenT< reco::SimToRecoCollection > simToRecoAssociationToken_
MonitorElement * meTrackEtaTot_
static constexpr double etaMatchCut_
static constexpr double zETL_
static constexpr double depositETLthreshold_
MonitorElement * meTrackMatchedTPEffEtaEtl2Mtd_
unsigned long long uint64_t
Definition: Time.h:13
MonitorElement * meMVATrackPullTot_
MonitorElement * meTrackMatchedTPEffPtMtd_
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
MonitorElement * meBTLTrackMatchedTPDPtvsPtGen_
void analyze(const edm::Event &, const edm::EventSetup &) override
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
double b
Definition: hdecay.h:120
edm::EDGetTokenT< CrossingFrame< PSimHit > > etlSimHitsToken_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
MonitorElement * meTrackSigmat0SafePid_
MonitorElement * meETLTrackEffEtaMtd_[2]
int zside() const
Definition: MTDDetId.h:61
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:664
MonitorElement * meTrackSigmat0Pid_
edm::EDGetTokenT< CrossingFrame< PSimHit > > btlSimHitsToken_
MonitorElement * meTrackMatchedTPmtdEffPtTot_
MonitorElement * meBTLTrackEffPtMtd_
MonitorElement * meTrackt0Src_
edm::EDGetTokenT< edm::ValueMap< float > > Sigmat0PidToken_
static constexpr double deltaDRcut_
MonitorElement * meTracktmtd_
MonitorElement * meBTLTrackEffPtTot_
MonitorElement * meTrackPathLenghtvsEta_
MonitorElement * meTrackMatchedTPEffPtTot_
Detector identifier class for the Endcap Timing Layer.
Definition: ETLDetId.h:16
HLT enums.
MonitorElement * meETLTrackMatchedTPPtResMtd_
int nDisc() const
Definition: ETLDetId.h:155
double a
Definition: hdecay.h:121
const bool mvaGenRecMatch(const HepMC::GenParticle &, const double &, const reco::TrackBase &, const bool &)
Detector identifier class for the Barrel Timing Layer. The crystal count must start from 0...
Definition: BTLDetId.h:19
MonitorElement * meETLTrackMatchedTPPtRatioGen_
MonitorElement * meExtraPtMtd_
ESTransientHandle< T > getTransientHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:141
MonitorElement * meExtraPhiAtBTLmatched_
edm::EDGetTokenT< reco::RecoToSimCollection > recoToSimAssociationToken_
Monte Carlo truth information used for tracking validation.
MonitorElement * meBTLTrackMatchedTPPtRatioGen_
edm::EDGetTokenT< edm::ValueMap< float > > Sigmat0SrcToken_
MonitorElement * meETLTrackEffPtMtd_[2]
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
MonitorElement * meExtraPtEtl2Mtd_
BTLDetId::CrysLayout crysLayoutFromTopoMode(const int &topoMode)
edm::ESGetToken< MTDGeometry, MTDDigiGeometryRecord > mtdgeoToken_
std::vector< TrackingParticle > TrackingParticleCollection
MonitorElement * meTrackPtTot_
MonitorElement * meMVATrackZposResTot_
Log< level::Warning, false > LogWarning
auto makeValid(const U &iOtherHandleType) noexcept(false)
Definition: ValidHandle.h:52
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
MonitorElement * meETLTrackMatchedTP2DPtvsPtMtd_
MonitorElement * meMVATrackMatchedEffPtMtd_
MonitorElement * meTrackSigmat0Src_
const std::vector< const DetLayer * > & allETLLayers() const
return the ETL DetLayers (endcap), -Z to +Z
Definition: Run.h:45
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magfieldToken_
edm::EDGetTokenT< edm::ValueMap< float > > Sigmat0SafePidToken_
MtdTracksValidation(const edm::ParameterSet &)
#define LogDebug(id)
edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > particleTableToken_
MonitorElement * meExtraPhiAtBTL_