CMS 3D CMS Logo

MtdTracksValidation.cc
Go to the documentation of this file.
1 #include <string>
2 
8 
11 
18 
24 
28 
36 
47 
58 
61 
68 
71 #include "HepMC/GenRanges.h"
72 #include "CLHEP/Units/PhysicalConstants.h"
73 #include "MTDHit.h"
74 
76 public:
77  explicit MtdTracksValidation(const edm::ParameterSet&);
78  ~MtdTracksValidation() override;
79 
80  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
81 
82 private:
83  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
84 
85  void analyze(const edm::Event&, const edm::EventSetup&) override;
86 
87  const std::pair<bool, bool> checkAcceptance(
88  const reco::Track&, const edm::Event&, const edm::EventSetup&, size_t&, float&, float&, float&, float&);
89 
90  const bool mvaGenSel(const HepMC::GenParticle&, const float&);
91  const bool mvaTPSel(const TrackingParticle&);
92  const bool mvaRecSel(const reco::TrackBase&, const reco::Vertex&, const double&, const double&);
93  const bool mvaGenRecMatch(const HepMC::GenParticle&, const double&, const reco::TrackBase&, const bool&);
95 
96  const unsigned long int uniqueId(const uint32_t x, const EncodedEventId& y) {
97  const uint64_t a = static_cast<uint64_t>(x);
98  const uint64_t b = static_cast<uint64_t>(y.rawId());
99 
100  if (x < y.rawId())
101  return (b << 32) | a;
102  else
103  return (a << 32) | b;
104  }
105 
106  bool isETL(const double eta) const { return (std::abs(eta) > trackMinEtlEta_) && (std::abs(eta) < trackMaxEtlEta_); }
107 
108  // ------------ member data ------------
109 
111  const float trackMinPt_;
112  const float trackMaxBtlEta_;
113  const float trackMinEtlEta_;
114  const float trackMaxEtlEta_;
115 
116  static constexpr double etacutGEN_ = 4.; // |eta| < 4;
117  static constexpr double etacutREC_ = 3.; // |eta| < 3;
118  static constexpr double pTcut_ = 0.7; // PT > 0.7 GeV
119  static constexpr double deltaZcut_ = 0.1; // dz separation 1 mm
120  static constexpr double deltaPTcut_ = 0.05; // dPT < 5%
121  static constexpr double deltaDRcut_ = 0.03; // DeltaR separation
122  static constexpr double depositBTLthreshold_ = 1; // threshold for energy deposit in BTL cell [MeV]
123  static constexpr double depositETLthreshold_ = 0.001; // threshold for energy deposit in ETL cell [MeV]
124  static constexpr double rBTL_ = 110.0;
125  static constexpr double zETL_ = 290.0;
126  static constexpr double etaMatchCut_ = 0.05;
127  static constexpr double cluDRradius_ = 0.05; // to cluster rechits around extrapolated track
128 
130 
133 
137 
146 
149 
159 
166 
175 
187 
199 
225 
229 };
230 
231 // ------------ constructor and destructor --------------
233  : folder_(iConfig.getParameter<std::string>("folder")),
234  trackMinPt_(iConfig.getParameter<double>("trackMinimumPt")),
235  trackMaxBtlEta_(iConfig.getParameter<double>("trackMaximumBtlEta")),
236  trackMinEtlEta_(iConfig.getParameter<double>("trackMinimumEtlEta")),
237  trackMaxEtlEta_(iConfig.getParameter<double>("trackMaximumEtlEta")),
238  optionalPlots_(iConfig.getUntrackedParameter<bool>("optionalPlots")) {
239  GenRecTrackToken_ = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("inputTagG"));
240  RecTrackToken_ = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("inputTagT"));
241  RecVertexToken_ = consumes<std::vector<reco::Vertex>>(iConfig.getParameter<edm::InputTag>("inputTagV"));
242  HepMCProductToken_ = consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("inputTagH"));
244  consumes<TrackingParticleCollection>(iConfig.getParameter<edm::InputTag>("SimTag"));
246  consumes<reco::SimToRecoCollection>(iConfig.getParameter<edm::InputTag>("TPtoRecoTrackAssoc"));
248  consumes<reco::RecoToSimCollection>(iConfig.getParameter<edm::InputTag>("TPtoRecoTrackAssoc"));
249  btlSimHitsToken_ = consumes<CrossingFrame<PSimHit>>(iConfig.getParameter<edm::InputTag>("btlSimHits"));
250  etlSimHitsToken_ = consumes<CrossingFrame<PSimHit>>(iConfig.getParameter<edm::InputTag>("etlSimHits"));
251  btlRecHitsToken_ = consumes<FTLRecHitCollection>(iConfig.getParameter<edm::InputTag>("btlRecHits"));
252  etlRecHitsToken_ = consumes<FTLRecHitCollection>(iConfig.getParameter<edm::InputTag>("etlRecHits"));
253  trackAssocToken_ = consumes<edm::ValueMap<int>>(iConfig.getParameter<edm::InputTag>("trackAssocSrc"));
254  pathLengthToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("pathLengthSrc"));
255  tmtdToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("tmtd"));
256  SigmatmtdToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmatmtd"));
257  t0SrcToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("t0Src"));
258  Sigmat0SrcToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmat0Src"));
259  t0PidToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("t0PID"));
260  Sigmat0PidToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmat0PID"));
261  t0SafePidToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("t0SafePID"));
262  Sigmat0SafePidToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmat0SafePID"));
263  trackMVAQualToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("trackMVAQual"));
264  mtdgeoToken_ = esConsumes<MTDGeometry, MTDDigiGeometryRecord>();
265  mtdtopoToken_ = esConsumes<MTDTopology, MTDTopologyRcd>();
266  mtdlayerToken_ = esConsumes<MTDDetLayerGeometry, MTDRecoGeometryRecord>();
267  magfieldToken_ = esConsumes<MagneticField, IdealMagneticFieldRecord>();
268  builderToken_ = esConsumes<TransientTrackBuilder, TransientTrackRecord>(edm::ESInputTag("", "TransientTrackBuilder"));
269  particleTableToken_ = esConsumes<HepPDT::ParticleDataTable, edm::DefaultRecord>();
270 }
271 
273 
274 // ------------ method called for each event ------------
276  using namespace edm;
277  using namespace geant_units::operators;
278  using namespace std;
279 
280  auto GenRecTrackHandle = makeValid(iEvent.getHandle(GenRecTrackToken_));
281  auto RecVertexHandle = makeValid(iEvent.getHandle(RecVertexToken_));
282 
283  std::unordered_map<uint32_t, MTDHit> m_btlHits;
284  std::unordered_map<uint32_t, MTDHit> m_etlHits;
285  std::unordered_map<uint32_t, std::set<unsigned long int>> m_btlTrkPerCell;
286  std::unordered_map<uint32_t, std::set<unsigned long int>> m_etlTrkPerCell;
287  std::map<TrackingParticleRef, std::vector<uint32_t>> m_tp2detid;
288 
289  const auto& tMtd = iEvent.get(tmtdToken_);
290  const auto& SigmatMtd = iEvent.get(SigmatmtdToken_);
291  const auto& t0Src = iEvent.get(t0SrcToken_);
292  const auto& Sigmat0Src = iEvent.get(Sigmat0SrcToken_);
293  const auto& t0Pid = iEvent.get(t0PidToken_);
294  const auto& Sigmat0Pid = iEvent.get(Sigmat0PidToken_);
295  const auto& t0Safe = iEvent.get(t0SafePidToken_);
296  const auto& Sigmat0Safe = iEvent.get(Sigmat0SafePidToken_);
297  const auto& mtdQualMVA = iEvent.get(trackMVAQualToken_);
298  const auto& trackAssoc = iEvent.get(trackAssocToken_);
299  const auto& pathLength = iEvent.get(pathLengthToken_);
300 
301  const auto& primRecoVtx = *(RecVertexHandle.product()->begin());
302 
303  // generator level information (HepMC format)
304  auto GenEventHandle = makeValid(iEvent.getHandle(HepMCProductToken_));
305  const HepMC::GenEvent* mc = GenEventHandle->GetEvent();
306  double zsim = convertMmToCm((*(mc->vertices_begin()))->position().z());
307  double tsim = (*(mc->vertices_begin()))->position().t() * CLHEP::mm / CLHEP::c_light;
308 
309  auto pdt = iSetup.getHandle(particleTableToken_);
310  const HepPDT::ParticleDataTable* pdTable = pdt.product();
311 
312  auto simToRecoH = makeValid(iEvent.getHandle(simToRecoAssociationToken_));
313  s2r_ = simToRecoH.product();
314 
315  auto recoToSimH = makeValid(iEvent.getHandle(recoToSimAssociationToken_));
316  r2s_ = recoToSimH.product();
317 
318  //Fill maps with simhits accumulated per DetId
319 
320  auto btlSimHitsHandle = makeValid(iEvent.getHandle(btlSimHitsToken_));
321  MixCollection<PSimHit> btlSimHits(btlSimHitsHandle.product());
322  for (auto const& simHit : btlSimHits) {
323  if (simHit.tof() < 0 || simHit.tof() > 25.)
324  continue;
325  DetId id = simHit.detUnitId();
326  auto const thisHId = uniqueId(simHit.trackId(), simHit.eventId());
327  m_btlTrkPerCell[id.rawId()].insert(thisHId);
328  auto simHitIt = m_btlHits.emplace(id.rawId(), MTDHit()).first;
329  // --- Accumulate the energy (in MeV) of SIM hits in the same detector cell
330  (simHitIt->second).energy += convertUnitsTo(0.001_MeV, simHit.energyLoss());
331  }
332 
333  auto etlSimHitsHandle = makeValid(iEvent.getHandle(etlSimHitsToken_));
334  MixCollection<PSimHit> etlSimHits(etlSimHitsHandle.product());
335  for (auto const& simHit : etlSimHits) {
336  if (simHit.tof() < 0 || simHit.tof() > 25.) {
337  continue;
338  }
339  DetId id = simHit.detUnitId();
340  auto const thisHId = uniqueId(simHit.trackId(), simHit.eventId());
341  m_etlTrkPerCell[id.rawId()].insert(thisHId);
342  auto simHitIt = m_etlHits.emplace(id.rawId(), MTDHit()).first;
343  // --- Accumulate the energy (in MeV) of SIM hits in the same detector cell
344  (simHitIt->second).energy += convertUnitsTo(0.001_MeV, simHit.energyLoss());
345  }
346 
347  //Fill map of DetId per ref to TP
348 
349  auto tpHandle = makeValid(iEvent.getHandle(trackingParticleCollectionToken_));
350  TrackingParticleCollection tpColl = *(tpHandle.product());
351  size_t tpindex(0);
352  for (auto tp = tpColl.begin(); tp != tpColl.end(); tp++, ++tpindex) {
354  if (tp->eventId().bunchCrossing() == 0 && tp->eventId().event() == 0) {
355  if (!mvaTPSel(*tp))
356  continue;
357  for (const auto& simTrk : tp->g4Tracks()) {
358  auto const thisTId = uniqueId(simTrk.trackId(), simTrk.eventId());
359  for (auto const& cell : m_btlTrkPerCell) {
360  if (m_btlHits[cell.first].energy < depositBTLthreshold_) {
361  continue;
362  }
363  for (auto const& simtrack : cell.second) {
364  if (thisTId == simtrack) {
365  m_tp2detid[tpref].emplace_back(cell.first);
366  break;
367  }
368  }
369  }
370  for (auto const& cell : m_etlTrkPerCell) {
371  if (m_etlHits[cell.first].energy < depositETLthreshold_) {
372  continue;
373  }
374  for (auto const& simtrack : cell.second) {
375  if (thisTId == simtrack) {
376  m_tp2detid[tpref].emplace_back(cell.first);
377  break;
378  }
379  }
380  }
381  }
382  }
383  }
384 
385  unsigned int index = 0;
386 
387  // flag to select events with reco vertex close to true simulated primary vertex, or PV fake (particle guns)
388  const bool isGoodVtx = std::abs(primRecoVtx.z() - zsim) < deltaZcut_ || primRecoVtx.isFake();
389 
390  // --- Loop over all RECO tracks ---
391  for (const auto& trackGen : *GenRecTrackHandle) {
392  const reco::TrackRef trackref(iEvent.getHandle(GenRecTrackToken_), index);
393  index++;
394 
395  if (trackAssoc[trackref] == -1) {
396  LogInfo("mtdTracks") << "Extended track not associated";
397  continue;
398  }
399 
400  const reco::TrackRef mtdTrackref = reco::TrackRef(iEvent.getHandle(RecTrackToken_), trackAssoc[trackref]);
401  const reco::Track& track = *mtdTrackref;
402 
403  bool isBTL = false;
404  bool isETL = false;
405  bool twoETLdiscs = false;
406  bool noCrack = std::abs(trackGen.eta()) < trackMaxBtlEta_ || std::abs(trackGen.eta()) > trackMinEtlEta_;
407 
408  if (track.pt() >= trackMinPt_ && std::abs(track.eta()) <= trackMaxEtlEta_) {
409  meTracktmtd_->Fill(tMtd[trackref]);
410  if (std::round(SigmatMtd[trackref] - Sigmat0Pid[trackref]) != 0) {
411  LogWarning("mtdTracks")
412  << "TimeError associated to refitted track is different from TimeError stored in tofPID "
413  "sigmat0 ValueMap: this should not happen";
414  }
415 
416  meTrackt0Src_->Fill(t0Src[trackref]);
417  meTrackSigmat0Src_->Fill(Sigmat0Src[trackref]);
418 
419  meTrackt0Pid_->Fill(t0Pid[trackref]);
420  meTrackSigmat0Pid_->Fill(Sigmat0Pid[trackref]);
421  meTrackt0SafePid_->Fill(t0Safe[trackref]);
422  meTrackSigmat0SafePid_->Fill(Sigmat0Safe[trackref]);
423  meTrackMVAQual_->Fill(mtdQualMVA[trackref]);
424 
425  meTrackPathLenghtvsEta_->Fill(std::abs(track.eta()), pathLength[trackref]);
426 
427  if (std::abs(track.eta()) < trackMaxBtlEta_) {
428  // --- all BTL tracks (with and without hit in MTD) ---
432 
433  bool MTDBtl = false;
434  int numMTDBtlvalidhits = 0;
435  for (const auto hit : track.recHits()) {
436  if (hit->isValid() == false)
437  continue;
438  MTDDetId Hit = hit->geographicalId();
439  if ((Hit.det() == 6) && (Hit.subdetId() == 1) && (Hit.mtdSubDetector() == 1)) {
440  MTDBtl = true;
441  numMTDBtlvalidhits++;
442  }
443  }
444  meTrackNumHits_->Fill(numMTDBtlvalidhits);
445  if (isBTL && Sigmat0Safe[trackref] < 0.) {
446  meTrackNumHitsNT_->Fill(numMTDBtlvalidhits);
447  }
448 
449  // --- keeping only tracks with last hit in MTD ---
450  if (MTDBtl == true) {
451  isBTL = true;
456  meBTLTrackPtRes_->Fill((trackGen.pt() - track.pt()) / trackGen.pt());
457  }
458  } //loop over (geometrical) BTL tracks
459 
460  else {
461  // --- all ETL tracks (with and without hit in MTD) ---
462  if ((track.eta() < -trackMinEtlEta_) && (track.eta() > -trackMaxEtlEta_)) {
463  meETLTrackEffEtaTot_[0]->Fill(track.eta());
464  meETLTrackEffPhiTot_[0]->Fill(track.phi());
465  meETLTrackEffPtTot_[0]->Fill(track.pt());
466  }
467 
468  if ((track.eta() > trackMinEtlEta_) && (track.eta() < trackMaxEtlEta_)) {
469  meETLTrackEffEtaTot_[1]->Fill(track.eta());
470  meETLTrackEffPhiTot_[1]->Fill(track.phi());
471  meETLTrackEffPtTot_[1]->Fill(track.pt());
472  }
473 
474  bool MTDEtlZnegD1 = false;
475  bool MTDEtlZnegD2 = false;
476  bool MTDEtlZposD1 = false;
477  bool MTDEtlZposD2 = false;
478  int numMTDEtlvalidhits = 0;
479  for (const auto hit : track.recHits()) {
480  if (hit->isValid() == false)
481  continue;
482  MTDDetId Hit = hit->geographicalId();
483  if ((Hit.det() == 6) && (Hit.subdetId() == 1) && (Hit.mtdSubDetector() == 2)) {
484  isETL = true;
485  ETLDetId ETLHit = hit->geographicalId();
486 
487  if ((ETLHit.zside() == -1) && (ETLHit.nDisc() == 1)) {
488  MTDEtlZnegD1 = true;
490  meETLTrackPtRes_->Fill((trackGen.pt() - track.pt()) / trackGen.pt());
491  numMTDEtlvalidhits++;
492  }
493  if ((ETLHit.zside() == -1) && (ETLHit.nDisc() == 2)) {
494  MTDEtlZnegD2 = true;
496  meETLTrackPtRes_->Fill((trackGen.pt() - track.pt()) / trackGen.pt());
497  numMTDEtlvalidhits++;
498  }
499  if ((ETLHit.zside() == 1) && (ETLHit.nDisc() == 1)) {
500  MTDEtlZposD1 = true;
502  meETLTrackPtRes_->Fill((trackGen.pt() - track.pt()) / trackGen.pt());
503  numMTDEtlvalidhits++;
504  }
505  if ((ETLHit.zside() == 1) && (ETLHit.nDisc() == 2)) {
506  MTDEtlZposD2 = true;
508  meETLTrackPtRes_->Fill((trackGen.pt() - track.pt()) / trackGen.pt());
509  numMTDEtlvalidhits++;
510  }
511  }
512  }
513  meTrackNumHits_->Fill(-numMTDEtlvalidhits);
514  if (isETL && Sigmat0Safe[trackref] < 0.) {
515  meTrackNumHitsNT_->Fill(-numMTDEtlvalidhits);
516  }
517 
518  // --- keeping only tracks with last hit in MTD ---
519  if ((track.eta() < -trackMinEtlEta_) && (track.eta() > -trackMaxEtlEta_)) {
520  twoETLdiscs = (MTDEtlZnegD1 == true) && (MTDEtlZnegD2 == true);
521  if ((MTDEtlZnegD1 == true) || (MTDEtlZnegD2 == true)) {
522  meETLTrackEffEtaMtd_[0]->Fill(track.eta());
523  meETLTrackEffPhiMtd_[0]->Fill(track.phi());
524  meETLTrackEffPtMtd_[0]->Fill(track.pt());
525  if (twoETLdiscs) {
526  meETLTrackEffEta2Mtd_[0]->Fill(track.eta());
527  meETLTrackEffPhi2Mtd_[0]->Fill(track.phi());
528  meETLTrackEffPt2Mtd_[0]->Fill(track.pt());
529  }
530  }
531  }
532  if ((track.eta() > trackMinEtlEta_) && (track.eta() < trackMaxEtlEta_)) {
533  twoETLdiscs = (MTDEtlZposD1 == true) && (MTDEtlZposD2 == true);
534  if ((MTDEtlZposD1 == true) || (MTDEtlZposD2 == true)) {
535  meETLTrackEffEtaMtd_[1]->Fill(track.eta());
536  meETLTrackEffPhiMtd_[1]->Fill(track.phi());
537  meETLTrackEffPtMtd_[1]->Fill(track.pt());
538  if (twoETLdiscs) {
539  meETLTrackEffEta2Mtd_[1]->Fill(track.eta());
540  meETLTrackEffPhi2Mtd_[1]->Fill(track.phi());
541  meETLTrackEffPt2Mtd_[1]->Fill(track.pt());
542  }
543  }
544  }
545  }
546 
547  // TrackingParticle based matching
548 
549  const reco::TrackBaseRef trkrefb(trackref);
550  auto tp_info = getMatchedTP(trkrefb);
551 
552  meTrackPtTot_->Fill(trackGen.pt());
553  meTrackEtaTot_->Fill(std::abs(trackGen.eta()));
554  if (tp_info != nullptr && mvaTPSel(**tp_info)) {
555  const bool withMTD = (m_tp2detid.find(*tp_info) != m_tp2detid.end());
556  if (noCrack) {
557  meTrackMatchedTPEffPtTot_->Fill(trackGen.pt());
558  if (withMTD) {
559  meTrackMatchedTPmtdEffPtTot_->Fill(trackGen.pt());
560  }
561  }
562  meTrackMatchedTPEffEtaTot_->Fill(std::abs(trackGen.eta()));
563  if (withMTD) {
564  meTrackMatchedTPmtdEffEtaTot_->Fill(std::abs(trackGen.eta()));
565  }
566  if (isBTL || isETL) {
567  if (noCrack) {
568  meTrackMatchedTPEffPtMtd_->Fill(trackGen.pt());
569  if (isBTL || twoETLdiscs) {
570  meTrackMatchedTPEffPtEtl2Mtd_->Fill(trackGen.pt());
571  }
572  if (withMTD) {
573  meTrackMatchedTPmtdEffPtMtd_->Fill(trackGen.pt());
574  }
575  }
576  meTrackMatchedTPEffEtaMtd_->Fill(std::abs(trackGen.eta()));
577  if (isBTL || twoETLdiscs) {
578  meTrackMatchedTPEffEtaEtl2Mtd_->Fill(std::abs(trackGen.eta()));
579  }
580  if (withMTD) {
581  meTrackMatchedTPmtdEffEtaMtd_->Fill(std::abs(trackGen.eta()));
582  }
583  }
584 
585  if (optionalPlots_) {
586  size_t nlayers(0);
587  float extrho(0.);
588  float exteta(0.);
589  float extphi(0.);
590  float selvar(0.);
591  auto accept = checkAcceptance(trackGen, iEvent, iSetup, nlayers, extrho, exteta, extphi, selvar);
592  if (accept.first && std::abs(exteta) < trackMaxBtlEta_) {
594  meExtraBTLeneInCone_->Fill(selvar);
595  }
596  if (accept.second) {
597  if (std::abs(exteta) < trackMaxBtlEta_) {
599  }
600  if (noCrack) {
601  meExtraPtMtd_->Fill(trackGen.pt());
602  if (nlayers == 2) {
603  meExtraPtEtl2Mtd_->Fill(trackGen.pt());
604  }
605  }
606  meExtraEtaMtd_->Fill(std::abs(trackGen.eta()));
607  if (nlayers == 2) {
608  meExtraEtaEtl2Mtd_->Fill(trackGen.eta());
609  }
610  }
611  }
612 
613  } // TP matching
614  }
615 
616  if (isGoodVtx) {
617  const bool vtxFake = primRecoVtx.isFake();
618 
619  if (mvaRecSel(trackGen, primRecoVtx, t0Safe[trackref], Sigmat0Safe[trackref])) {
620  // reco-gen matching used for MVA quality flag
621 
622  if (noCrack) {
623  meMVATrackEffPtTot_->Fill(trackGen.pt());
624  }
625  meMVATrackEffEtaTot_->Fill(std::abs(trackGen.eta()));
626 
627  double dZ = trackGen.vz() - zsim;
628  double dT(-9999.);
629  double pullT(-9999.);
630  if (Sigmat0Safe[trackref] != -1.) {
631  dT = t0Safe[trackref] - tsim;
632  pullT = dT / Sigmat0Safe[trackref];
633  }
634  for (const auto& genP : mc->particle_range()) {
635  // select status 1 genParticles and match them to the reconstructed track
636 
637  float charge = pdTable->particle(HepPDT::ParticleID(genP->pdg_id())) != nullptr
638  ? pdTable->particle(HepPDT::ParticleID(genP->pdg_id()))->charge()
639  : 0.f;
640  if (mvaGenSel(*genP, charge)) {
641  if (mvaGenRecMatch(*genP, zsim, trackGen, vtxFake)) {
643  if (noCrack) {
644  meMVATrackMatchedEffPtTot_->Fill(trackGen.pt());
645  }
646  meMVATrackMatchedEffEtaTot_->Fill(std::abs(trackGen.eta()));
647  if (isBTL || isETL) {
648  meMVATrackResTot_->Fill(dT);
649  meMVATrackPullTot_->Fill(pullT);
650  if (noCrack) {
651  meMVATrackMatchedEffPtMtd_->Fill(trackGen.pt());
652  }
653  meMVATrackMatchedEffEtaMtd_->Fill(std::abs(trackGen.eta()));
654  }
655  break;
656  }
657  }
658  }
659  }
660  } // MC truth matich analysis for good PV
661  } //RECO tracks loop
662 }
663 
664 const std::pair<bool, bool> MtdTracksValidation::checkAcceptance(const reco::Track& track,
665  const edm::Event& iEvent,
666  edm::EventSetup const& iSetup,
667  size_t& nlayers,
668  float& extrho,
669  float& exteta,
670  float& extphi,
671  float& selvar) {
672  bool isMatched(false);
673  nlayers = 0;
674  extrho = 0.;
675  exteta = -999.;
676  extphi = -999.;
677  selvar = 0.;
678 
679  auto geometryHandle = iSetup.getTransientHandle(mtdgeoToken_);
680  const MTDGeometry* geom = geometryHandle.product();
681  auto topologyHandle = iSetup.getTransientHandle(mtdtopoToken_);
682  const MTDTopology* topology = topologyHandle.product();
683 
684  auto layerHandle = iSetup.getTransientHandle(mtdlayerToken_);
685  const MTDDetLayerGeometry* layerGeo = layerHandle.product();
686 
687  auto magfieldHandle = iSetup.getTransientHandle(magfieldToken_);
688  const MagneticField* mfield = magfieldHandle.product();
689 
690  auto ttrackBuilder = iSetup.getTransientHandle(builderToken_);
691 
692  auto tTrack = ttrackBuilder->build(track);
693  TrajectoryStateOnSurface tsos = tTrack.outermostMeasurementState();
694  float theMaxChi2 = 500.;
695  float theNSigma = 10.;
696  std::unique_ptr<MeasurementEstimator> theEstimator =
697  std::make_unique<Chi2MeasurementEstimator>(theMaxChi2, theNSigma);
699 
700  auto btlRecHitsHandle = makeValid(iEvent.getHandle(btlRecHitsToken_));
701  auto etlRecHitsHandle = makeValid(iEvent.getHandle(etlRecHitsToken_));
702 
703  edm::LogVerbatim("MtdTracksValidation") << "New Track, pt= " << track.pt() << " eta= " << track.eta();
704 
705  //try BTL
706  bool inBTL = false;
707  float eneSum(0.);
708  const std::vector<const DetLayer*>& layersBTL = layerGeo->allBTLLayers();
709  for (const DetLayer* ilay : layersBTL) {
710  std::pair<bool, TrajectoryStateOnSurface> comp = ilay->compatible(tsos, prop, *theEstimator);
711  if (!comp.first)
712  continue;
713  if (!inBTL) {
714  inBTL = true;
715  extrho = comp.second.globalPosition().perp();
716  exteta = comp.second.globalPosition().eta();
717  extphi = comp.second.globalPosition().phi();
718  edm::LogVerbatim("MtdTracksValidation")
719  << "Extrapolation at BTL surface, rho= " << extrho << " eta= " << exteta << " phi= " << extphi;
720  }
721  std::vector<DetLayer::DetWithState> compDets = ilay->compatibleDets(tsos, prop, *theEstimator);
722  for (const auto& detWithState : compDets) {
723  const auto& det = detWithState.first;
724 
725  // loop on compatible rechits and check energy in a fixed size cone around the extrapolation point
726 
727  edm::LogVerbatim("MtdTracksValidation")
728  << "DetId= " << det->geographicalId().rawId() << " gp= " << detWithState.second.globalPosition().x() << " "
729  << detWithState.second.globalPosition().y() << " " << detWithState.second.globalPosition().z()
730  << " rho= " << detWithState.second.globalPosition().perp()
731  << " eta= " << detWithState.second.globalPosition().eta()
732  << " phi= " << detWithState.second.globalPosition().phi();
733 
734  for (const auto& recHit : *btlRecHitsHandle) {
735  BTLDetId detId = recHit.id();
736  DetId geoId = detId.geographicalId(MTDTopologyMode::crysLayoutFromTopoMode(topology->getMTDTopologyMode()));
737  const MTDGeomDet* thedet = geom->idToDet(geoId);
738  if (thedet == nullptr)
739  throw cms::Exception("MtdTracksValidation") << "GeographicalID: " << std::hex << geoId.rawId() << " ("
740  << detId.rawId() << ") is invalid!" << std::dec << std::endl;
741  if (geoId == det->geographicalId()) {
742  const ProxyMTDTopology& topoproxy = static_cast<const ProxyMTDTopology&>(thedet->topology());
743  const RectangularMTDTopology& topo = static_cast<const RectangularMTDTopology&>(topoproxy.specificTopology());
744 
745  Local3DPoint local_point(0., 0., 0.);
746  local_point = topo.pixelToModuleLocalPoint(local_point, detId.row(topo.nrows()), detId.column(topo.nrows()));
747  const auto& global_point = thedet->toGlobal(local_point);
748  edm::LogVerbatim("MtdTracksValidation")
749  << "Hit id= " << detId.rawId() << " ene= " << recHit.energy()
750  << " dr= " << reco::deltaR(global_point, detWithState.second.globalPosition());
751  if (reco::deltaR(global_point, detWithState.second.globalPosition()) < cluDRradius_) {
752  eneSum += recHit.energy();
753  //extrho = detWithState.second.globalPosition().perp();
754  //exteta = detWithState.second.globalPosition().eta();
755  //extphi = detWithState.second.globalPosition().phi();
756  }
757  }
758  }
759  }
760  if (eneSum > depositBTLthreshold_) {
761  nlayers++;
762  selvar = eneSum;
763  isMatched = true;
764  edm::LogVerbatim("MtdTracksValidation") << "BTL matched, energy= " << eneSum << " #layers= " << nlayers;
765  }
766  }
767  if (inBTL) {
768  return std::make_pair(inBTL, isMatched);
769  }
770 
771  //try ETL
772  bool inETL = false;
773  const std::vector<const DetLayer*>& layersETL = layerGeo->allETLLayers();
774  for (const DetLayer* ilay : layersETL) {
775  size_t hcount(0);
776  const BoundDisk& disk = static_cast<const MTDRingForwardDoubleLayer*>(ilay)->specificSurface();
777  const double diskZ = disk.position().z();
778  if (tsos.globalPosition().z() * diskZ < 0)
779  continue; // only propagate to the disk that's on the same side
780  std::pair<bool, TrajectoryStateOnSurface> comp = ilay->compatible(tsos, prop, *theEstimator);
781  if (!comp.first)
782  continue;
783  if (!inETL) {
784  inETL = true;
785  extrho = comp.second.globalPosition().perp();
786  exteta = comp.second.globalPosition().eta();
787  extphi = comp.second.globalPosition().phi();
788  }
789  edm::LogVerbatim("MtdTracksValidation")
790  << "Extrapolation at ETL surface, rho= " << extrho << " eta= " << exteta << " phi= " << extphi;
791  std::vector<DetLayer::DetWithState> compDets = ilay->compatibleDets(tsos, prop, *theEstimator);
792  for (const auto& detWithState : compDets) {
793  const auto& det = detWithState.first;
794 
795  // loop on compatible rechits and check hits in a fixed size cone around the extrapolation point
796 
797  edm::LogVerbatim("MtdTracksValidation")
798  << "DetId= " << det->geographicalId().rawId() << " gp= " << detWithState.second.globalPosition().x() << " "
799  << detWithState.second.globalPosition().y() << " " << detWithState.second.globalPosition().z()
800  << " rho= " << detWithState.second.globalPosition().perp()
801  << " eta= " << detWithState.second.globalPosition().eta()
802  << " phi= " << detWithState.second.globalPosition().phi();
803 
804  for (const auto& recHit : *etlRecHitsHandle) {
805  ETLDetId detId = recHit.id();
806  DetId geoId = detId.geographicalId();
807  const MTDGeomDet* thedet = geom->idToDet(geoId);
808  if (thedet == nullptr)
809  throw cms::Exception("MtdTracksValidation") << "GeographicalID: " << std::hex << geoId.rawId() << " ("
810  << detId.rawId() << ") is invalid!" << std::dec << std::endl;
811  if (geoId == det->geographicalId()) {
812  const ProxyMTDTopology& topoproxy = static_cast<const ProxyMTDTopology&>(thedet->topology());
813  const RectangularMTDTopology& topo = static_cast<const RectangularMTDTopology&>(topoproxy.specificTopology());
814 
815  Local3DPoint local_point(topo.localX(recHit.row()), topo.localY(recHit.column()), 0.);
816  const auto& global_point = thedet->toGlobal(local_point);
817  edm::LogVerbatim("MtdTracksValidation")
818  << "Hit id= " << detId.rawId() << " time= " << recHit.time()
819  << " dr= " << reco::deltaR(global_point, detWithState.second.globalPosition());
820  if (reco::deltaR(global_point, detWithState.second.globalPosition()) < cluDRradius_) {
821  hcount++;
822  if (hcount == 1) {
823  //extrho = detWithState.second.globalPosition().perp();
824  //exteta = detWithState.second.globalPosition().eta();
825  //extphi = detWithState.second.globalPosition().phi();
826  }
827  }
828  }
829  }
830  }
831  if (hcount > 0) {
832  nlayers++;
833  selvar = (float)hcount;
834  isMatched = true;
835  edm::LogVerbatim("MtdTracksValidation") << "ETL matched, counts= " << hcount << " #layers= " << nlayers;
836  }
837  }
838 
839  if (!inBTL && !inETL) {
840  edm::LogVerbatim("MtdTracksValidation")
841  << "Track not extrapolating to MTD: pt= " << track.pt() << " eta= " << track.eta() << " phi= " << track.phi()
842  << " vz= " << track.vz() << " vxy= " << std::sqrt(track.vx() * track.vx() + track.vy() * track.vy());
843  }
844  return std::make_pair(inETL, isMatched);
845 }
846 
847 // ------------ method for histogram booking ------------
849  ibook.setCurrentFolder(folder_);
850 
851  // histogram booking
852  meBTLTrackRPTime_ = ibook.book1D("TrackBTLRPTime", "Track t0 with respect to R.P.;t0 [ns]", 100, -1, 3);
853  meBTLTrackEffEtaTot_ = ibook.book1D("TrackBTLEffEtaTot", "Track efficiency vs eta (Tot);#eta_{RECO}", 100, -1.6, 1.6);
855  ibook.book1D("TrackBTLEffPhiTot", "Track efficiency vs phi (Tot);#phi_{RECO} [rad]", 100, -3.2, 3.2);
856  meBTLTrackEffPtTot_ = ibook.book1D("TrackBTLEffPtTot", "Track efficiency vs pt (Tot);pt_{RECO} [GeV]", 50, 0, 10);
857  meBTLTrackEffEtaMtd_ = ibook.book1D("TrackBTLEffEtaMtd", "Track efficiency vs eta (Mtd);#eta_{RECO}", 100, -1.6, 1.6);
859  ibook.book1D("TrackBTLEffPhiMtd", "Track efficiency vs phi (Mtd);#phi_{RECO} [rad]", 100, -3.2, 3.2);
860  meBTLTrackEffPtMtd_ = ibook.book1D("TrackBTLEffPtMtd", "Track efficiency vs pt (Mtd);pt_{RECO} [GeV]", 50, 0, 10);
862  ibook.book1D("TrackBTLPtRes", "Track pT resolution ;pT_{Gentrack}-pT_{MTDtrack}/pT_{Gentrack} ", 100, -0.1, 0.1);
863  meETLTrackRPTime_ = ibook.book1D("TrackETLRPTime", "Track t0 with respect to R.P.;t0 [ns]", 100, -1, 3);
865  ibook.book1D("TrackETLEffEtaTotZneg", "Track efficiency vs eta (Tot) (-Z);#eta_{RECO}", 100, -3.2, -1.4);
867  ibook.book1D("TrackETLEffEtaTotZpos", "Track efficiency vs eta (Tot) (+Z);#eta_{RECO}", 100, 1.4, 3.2);
869  ibook.book1D("TrackETLEffPhiTotZneg", "Track efficiency vs phi (Tot) (-Z);#phi_{RECO} [rad]", 100, -3.2, 3.2);
871  ibook.book1D("TrackETLEffPhiTotZpos", "Track efficiency vs phi (Tot) (+Z);#phi_{RECO} [rad]", 100, -3.2, 3.2);
873  ibook.book1D("TrackETLEffPtTotZneg", "Track efficiency vs pt (Tot) (-Z);pt_{RECO} [GeV]", 50, 0, 10);
875  ibook.book1D("TrackETLEffPtTotZpos", "Track efficiency vs pt (Tot) (+Z);pt_{RECO} [GeV]", 50, 0, 10);
877  ibook.book1D("TrackETLEffEtaMtdZneg", "Track efficiency vs eta (Mtd) (-Z);#eta_{RECO}", 100, -3.2, -1.4);
879  ibook.book1D("TrackETLEffEtaMtdZpos", "Track efficiency vs eta (Mtd) (+Z);#eta_{RECO}", 100, 1.4, 3.2);
881  ibook.book1D("TrackETLEffPhiMtdZneg", "Track efficiency vs phi (Mtd) (-Z);#phi_{RECO} [rad]", 100, -3.2, 3.2);
883  ibook.book1D("TrackETLEffPhiMtdZpos", "Track efficiency vs phi (Mtd) (+Z);#phi_{RECO} [rad]", 100, -3.2, 3.2);
885  ibook.book1D("TrackETLEffPtMtdZneg", "Track efficiency vs pt (Mtd) (-Z);pt_{RECO} [GeV]", 50, 0, 10);
887  ibook.book1D("TrackETLEffPtMtdZpos", "Track efficiency vs pt (Mtd) (+Z);pt_{RECO} [GeV]", 50, 0, 10);
889  ibook.book1D("TrackETLEffEta2MtdZneg", "Track efficiency vs eta (Mtd 2 hit) (-Z);#eta_{RECO}", 100, -3.2, -1.4);
891  ibook.book1D("TrackETLEffEta2MtdZpos", "Track efficiency vs eta (Mtd 2 hit) (+Z);#eta_{RECO}", 100, 1.4, 3.2);
892  meETLTrackEffPhi2Mtd_[0] = ibook.book1D(
893  "TrackETLEffPhi2MtdZneg", "Track efficiency vs phi (Mtd 2 hit) (-Z);#phi_{RECO} [rad]", 100, -3.2, 3.2);
894  meETLTrackEffPhi2Mtd_[1] = ibook.book1D(
895  "TrackETLEffPhi2MtdZpos", "Track efficiency vs phi (Mtd 2 hit) (+Z);#phi_{RECO} [rad]", 100, -3.2, 3.2);
897  ibook.book1D("TrackETLEffPt2MtdZneg", "Track efficiency vs pt (Mtd 2 hit) (-Z);pt_{RECO} [GeV]", 50, 0, 10);
899  ibook.book1D("TrackETLEffPt2MtdZpos", "Track efficiency vs pt (Mtd 2 hit) (+Z);pt_{RECO} [GeV]", 50, 0, 10);
901  ibook.book1D("TrackETLPtRes", "Track pT resolution;pT_{Gentrack}-pT_{MTDtrack}/pT_{Gentrack} ", 100, -0.1, 0.1);
902 
903  meTracktmtd_ = ibook.book1D("Tracktmtd", "Track time from TrackExtenderWithMTD;tmtd [ns]", 150, 1, 16);
904  meTrackt0Src_ = ibook.book1D("Trackt0Src", "Track time from TrackExtenderWithMTD;t0Src [ns]", 100, -1.5, 1.5);
906  ibook.book1D("TrackSigmat0Src", "Time Error from TrackExtenderWithMTD; #sigma_{t0Src} [ns]", 100, 0, 0.1);
907 
908  meTrackt0Pid_ = ibook.book1D("Trackt0Pid", "Track t0 as stored in TofPid;t0 [ns]", 100, -1, 1);
909  meTrackSigmat0Pid_ = ibook.book1D("TrackSigmat0Pid", "Sigmat0 as stored in TofPid; #sigma_{t0} [ns]", 100, 0, 0.1);
910  meTrackt0SafePid_ = ibook.book1D("Trackt0SafePID", "Track t0 Safe as stored in TofPid;t0 [ns]", 100, -1, 1);
912  ibook.book1D("TrackSigmat0SafePID", "Sigmat0 Safe as stored in TofPid; #sigma_{t0} [ns]", 100, 0, 0.1);
913  meTrackNumHits_ = ibook.book1D("TrackNumHits", "Number of valid MTD hits per track ; Number of hits", 10, -5, 5);
914  meTrackNumHitsNT_ = ibook.book1D(
915  "TrackNumHitsNT", "Number of valid MTD hits per track no time associated; Number of hits", 10, -5, 5);
916  meTrackMVAQual_ = ibook.book1D("TrackMVAQual", "Track MVA Quality as stored in Value Map ; MVAQual", 100, 0, 1);
918  "TrackPathLenghtvsEta", "MTD Track pathlength vs MTD track Eta;|#eta|;Pathlength", 100, 0, 3.2, 100.0, 400.0, "S");
919 
920  meMVATrackEffPtTot_ = ibook.book1D("MVAEffPtTot", "Pt of tracks associated to LV; track pt [GeV] ", 110, 0., 11.);
922  ibook.book1D("MVAMatchedEffPtTot", "Pt of tracks associated to LV matched to GEN; track pt [GeV] ", 110, 0., 11.);
924  "MVAMatchedEffPtMtd", "Pt of tracks associated to LV matched to GEN with time; track pt [GeV] ", 110, 0., 11.);
925 
926  if (optionalPlots_) {
927  meExtraPtMtd_ = ibook.book1D("ExtraPtMtd", "Pt of tracks extrapolated to hits; track pt [GeV] ", 110, 0., 11.);
928  meExtraPtEtl2Mtd_ = ibook.book1D(
929  "ExtraPtEtl2Mtd", "Pt of tracks extrapolated to hits, 2 ETL layers; track pt [GeV] ", 110, 0., 11.);
930  }
931 
932  meTrackPtTot_ = ibook.book1D("TrackPtTot", "Pt of tracks ; track pt [GeV] ", 110, 0., 11.);
934  ibook.book1D("MatchedTPEffPtTot", "Pt of tracks matched to TP; track pt [GeV] ", 110, 0., 11.);
936  ibook.book1D("MatchedTPEffPtMtd", "Pt of tracks matched to TP with time; track pt [GeV] ", 110, 0., 11.);
938  "MatchedTPEffPtEtl2Mtd", "Pt of tracks matched to TP with time, 2 ETL hits; track pt [GeV] ", 110, 0., 11.);
939 
941  ibook.book1D("MatchedTPmtdEffPtTot", "Pt of tracks matched to TP-mtd hit; track pt [GeV] ", 110, 0., 11.);
943  "MatchedTPmtdEffPtMtd", "Pt of tracks matched to TP-mtd hit with time; track pt [GeV] ", 110, 0., 11.);
944 
945  meMVATrackEffEtaTot_ = ibook.book1D("MVAEffEtaTot", "Eta of tracks associated to LV; track eta ", 66, 0., 3.3);
947  ibook.book1D("MVAMatchedEffEtaTot", "Eta of tracks associated to LV matched to GEN; track eta ", 66, 0., 3.3);
949  "MVAMatchedEffEtaMtd", "Eta of tracks associated to LV matched to GEN with time; track eta ", 66, 0., 3.3);
950 
951  if (optionalPlots_) {
952  meExtraEtaMtd_ = ibook.book1D("ExtraEtaMtd", "Eta of tracks extrapolated to hits; track eta ", 66, 0., 3.3);
954  ibook.book1D("ExtraEtaEtl2Mtd", "Eta of tracks extrapolated to hits, 2 ETL layers; track eta ", 66, 0., 3.3);
955  }
956 
957  meTrackEtaTot_ = ibook.book1D("TrackEtaTot", "Eta of tracks ; track eta ", 66, 0., 3.3);
959  ibook.book1D("MatchedTPEffEtaTot", "Eta of tracks matched to TP; track eta ", 66, 0., 3.3);
960  meMVATrackEffEtaTot_ = ibook.book1D("MVAEffEtaTot", "Eta of tracks ; track eta ", 66, 0., 3.3);
962  ibook.book1D("MatchedTPEffEtaMtd", "Eta of tracks matched to TP with time; track eta ", 66, 0., 3.3);
964  "MatchedTPEffEtaEtl2Mtd", "Eta of tracks matched to TP with time, 2 ETL hits; track eta ", 66, 0., 3.3);
965 
967  ibook.book1D("MatchedTPmtdEffEtaTot", "Eta of tracks matched to TP-mtd hit; track eta ", 66, 0., 3.3);
969  ibook.book1D("MatchedTPmtdEffEtaMtd", "Eta of tracks matched to TP-mtd hit with time; track eta ", 66, 0., 3.3);
970 
971  meMVATrackResTot_ = ibook.book1D(
972  "MVATrackRes", "t_{rec} - t_{sim} for LV associated tracks; t_{rec} - t_{sim} [ns] ", 120, -0.15, 0.15);
974  ibook.book1D("MVATrackPull", "Pull for associated tracks; (t_{rec}-t_{sim})/#sigma_{t}", 50, -5., 5.);
976  "MVATrackZposResTot", "Z_{PCA} - Z_{sim} for associated tracks;Z_{PCA} - Z_{sim} [cm] ", 100, -0.1, 0.1);
977 
978  if (optionalPlots_) {
980  ibook.book1D("ExtraPhiAtBTL", "Phi at BTL surface of extrapolated tracks; phi [deg]", 720, -180., 180.);
982  ibook.book1D("ExtraPhiAtBTLmatched",
983  "Phi at BTL surface of extrapolated tracksi matched with BTL hits; phi [deg]",
984  720,
985  -180.,
986  180.);
988  "ExtraBTLeneInCone", "BTL reconstructed energy in cone arounnd extrapolated track; E [MeV]", 100, 0., 50.);
989  }
990 }
991 
992 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
993 
996 
997  desc.add<std::string>("folder", "MTD/Tracks");
998  desc.add<edm::InputTag>("inputTagG", edm::InputTag("generalTracks"));
999  desc.add<edm::InputTag>("inputTagT", edm::InputTag("trackExtenderWithMTD"));
1000  desc.add<edm::InputTag>("inputTagV", edm::InputTag("offlinePrimaryVertices4D"));
1001  desc.add<edm::InputTag>("inputTagH", edm::InputTag("generatorSmeared"));
1002  desc.add<edm::InputTag>("SimTag", edm::InputTag("mix", "MergedTrackTruth"));
1003  desc.add<edm::InputTag>("TPtoRecoTrackAssoc", edm::InputTag("trackingParticleRecoTrackAsssociation"));
1004  desc.add<edm::InputTag>("btlSimHits", edm::InputTag("mix", "g4SimHitsFastTimerHitsBarrel"));
1005  desc.add<edm::InputTag>("etlSimHits", edm::InputTag("mix", "g4SimHitsFastTimerHitsEndcap"));
1006  desc.add<edm::InputTag>("btlRecHits", edm::InputTag("mtdRecHits", "FTLBarrel"));
1007  desc.add<edm::InputTag>("etlRecHits", edm::InputTag("mtdRecHits", "FTLEndcap"));
1008  desc.add<edm::InputTag>("tmtd", edm::InputTag("trackExtenderWithMTD:generalTracktmtd"));
1009  desc.add<edm::InputTag>("sigmatmtd", edm::InputTag("trackExtenderWithMTD:generalTracksigmatmtd"));
1010  desc.add<edm::InputTag>("t0Src", edm::InputTag("trackExtenderWithMTD:generalTrackt0"));
1011  desc.add<edm::InputTag>("sigmat0Src", edm::InputTag("trackExtenderWithMTD:generalTracksigmat0"));
1012  desc.add<edm::InputTag>("trackAssocSrc", edm::InputTag("trackExtenderWithMTD:generalTrackassoc"))
1013  ->setComment("Association between General and MTD Extended tracks");
1014  desc.add<edm::InputTag>("pathLengthSrc", edm::InputTag("trackExtenderWithMTD:generalTrackPathLength"));
1015  desc.add<edm::InputTag>("t0SafePID", edm::InputTag("tofPID:t0safe"));
1016  desc.add<edm::InputTag>("sigmat0SafePID", edm::InputTag("tofPID:sigmat0safe"));
1017  desc.add<edm::InputTag>("sigmat0PID", edm::InputTag("tofPID:sigmat0"));
1018  desc.add<edm::InputTag>("t0PID", edm::InputTag("tofPID:t0"));
1019  desc.add<edm::InputTag>("trackMVAQual", edm::InputTag("mtdTrackQualityMVA:mtdQualMVA"));
1020  desc.add<double>("trackMinimumPt", 0.7); // [GeV]
1021  desc.add<double>("trackMaximumBtlEta", 1.5);
1022  desc.add<double>("trackMinimumEtlEta", 1.6);
1023  desc.add<double>("trackMaximumEtlEta", 3.);
1024  desc.addUntracked<bool>("optionalPlots", true);
1025 
1026  descriptions.add("mtdTracksValid", desc);
1027 }
1028 
1030  bool match = false;
1031  if (gp.status() != 1) {
1032  return match;
1033  }
1034  match = charge != 0.f && gp.momentum().perp() > pTcut_ && std::abs(gp.momentum().eta()) < etacutGEN_;
1035  return match;
1036 }
1037 
1039  bool match = false;
1040  if (tp.status() != 1) {
1041  return match;
1042  }
1043  auto x_pv = tp.parentVertex()->position().x();
1044  auto y_pv = tp.parentVertex()->position().y();
1045  auto z_pv = tp.parentVertex()->position().z();
1046 
1047  auto r_pv = std::sqrt(x_pv * x_pv + y_pv * y_pv);
1048 
1049  match = tp.charge() != 0 && tp.pt() > pTcut_ && std::abs(tp.eta()) < etacutGEN_ && r_pv < rBTL_ && z_pv < zETL_;
1050  return match;
1051 }
1052 
1054  const reco::Vertex& vtx,
1055  const double& t0,
1056  const double& st0) {
1057  bool match = false;
1058  match = trk.pt() > pTcut_ && std::abs(trk.eta()) < etacutREC_ &&
1059  (std::abs(trk.vz() - vtx.z()) <= deltaZcut_ || vtx.isFake());
1060  if (st0 > 0.) {
1061  match = match && std::abs(t0 - vtx.t()) < 3. * st0;
1062  }
1063  return match;
1064 }
1065 
1067  const double& zsim,
1068  const reco::TrackBase& trk,
1069  const bool& vtxFake) {
1070  bool match = false;
1071  double dR = reco::deltaR(genP.momentum(), trk.momentum());
1072  double genPT = genP.momentum().perp();
1073  match = std::abs(genPT - trk.pt()) < trk.pt() * deltaPTcut_ && dR < deltaDRcut_ &&
1074  (std::abs(trk.vz() - zsim) < deltaZcut_ || vtxFake);
1075  return match;
1076 }
1077 
1079  auto found = r2s_->find(recoTrack);
1080 
1081  // reco track not matched to any TP
1082  if (found == r2s_->end())
1083  return nullptr;
1084 
1085  //matched TP equal to any TP associated to in time events
1086  for (const auto& tp : found->val) {
1087  if (tp.first->eventId().bunchCrossing() == 0)
1088  return &tp.first;
1089  }
1090 
1091  // reco track not matched to any TP from vertex
1092  return nullptr;
1093 }
1094 
MonitorElement * meETLTrackEffEtaTot_[2]
edm::EDGetTokenT< edm::ValueMap< int > > trackAssocToken_
static constexpr double etacutGEN_
Log< level::Info, true > LogVerbatim
MonitorElement * meTrackMatchedTPEffEtaMtd_
MonitorElement * meBTLTrackEffPhiMtd_
int getMTDTopologyMode() const
Definition: MTDTopology.h:27
edm::EDGetTokenT< edm::ValueMap< float > > trackMVAQualToken_
edm::EDGetTokenT< edm::ValueMap< float > > pathLengthToken_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
MonitorElement * meMVATrackResTot_
MonitorElement * meTrackNumHits_
MonitorElement * meTrackt0Pid_
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 * meETLTrackRPTime_
edm::EDGetTokenT< edm::ValueMap< float > > t0SafePidToken_
std::string folder_
edm::EDGetTokenT< edm::ValueMap< float > > t0SrcToken_
static constexpr double rBTL_
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
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 * 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)
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 &)
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]
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:399
MonitorElement * meTrackt0SafePid_
T sqrt(T t)
Definition: SSEVec.h:19
MonitorElement * meBTLTrackEffEtaTot_
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)
#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 * 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
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:15
HLT enums.
int nDisc() const
Definition: ETLDetId.h:124
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 * 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.
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 * 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 &)
edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > particleTableToken_
MonitorElement * meExtraPhiAtBTL_