CMS 3D CMS Logo

MtdTracksValidation.cc
Go to the documentation of this file.
1 #include <string>
2 
8 
11 
17 
23 
27 
35 
42 
45 #include "HepMC/GenRanges.h"
46 #include "CLHEP/Units/PhysicalConstants.h"
47 
48 namespace {
49 
50  struct MTDHit {
51  float energy;
52  float time;
53  float x;
54  float y;
55  float z;
56  };
57 
58 } // namespace
59 
61 public:
62  explicit MtdTracksValidation(const edm::ParameterSet&);
63  ~MtdTracksValidation() override;
64 
65  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
66 
67 private:
68  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
69 
70  void analyze(const edm::Event&, const edm::EventSetup&) override;
71 
72  const bool mvaGenSel(const HepMC::GenParticle&, const float&);
73  const bool mvaTPSel(const TrackingParticle&);
74  const bool mvaRecSel(const reco::TrackBase&, const reco::Vertex&, const double&, const double&);
75  const bool mvaGenRecMatch(const HepMC::GenParticle&, const double&, const reco::TrackBase&, const bool&);
77  const bool tpWithMTD(const TrackingParticle&, const std::unordered_set<unsigned long int>&);
78 
79  const unsigned long int uniqueId(const uint32_t x, const EncodedEventId& y) {
80  const uint64_t a = static_cast<uint64_t>(x);
81  const uint64_t b = static_cast<uint64_t>(y.rawId());
82 
83  if (x < y.rawId())
84  return (b << 32) | a;
85  else
86  return (a << 32) | b;
87  }
88 
89  bool isETL(const double eta) const { return (std::abs(eta) > trackMinEtlEta_) && (std::abs(eta) < trackMaxEtlEta_); }
90 
91  // ------------ member data ------------
92 
94  const float trackMinPt_;
95  const float trackMaxBtlEta_;
96  const float trackMinEtlEta_;
97  const float trackMaxEtlEta_;
98 
99  static constexpr double etacutGEN_ = 4.; // |eta| < 4;
100  static constexpr double etacutREC_ = 3.; // |eta| < 3;
101  static constexpr double pTcut_ = 0.7; // PT > 0.7 GeV
102  static constexpr double deltaZcut_ = 0.1; // dz separation 1 mm
103  static constexpr double deltaPTcut_ = 0.05; // dPT < 5%
104  static constexpr double deltaDRcut_ = 0.03; // DeltaR separation
105  static constexpr double depositBTLthreshold_ = 1; // threshold for energy deposit in BTL cell [MeV]
106  static constexpr double depositETLthreshold_ = 0.001; // threshold for energy deposit in ETL cell [MeV]
107  static constexpr double rBTL_ = 110.0;
108  static constexpr double zETL_ = 290.0;
109  static constexpr double etaMatchCut_ = 0.05;
110 
112 
115 
119 
126 
129 
139 
143 
152 
164 
175 
195 
201 };
202 
203 // ------------ constructor and destructor --------------
205  : folder_(iConfig.getParameter<std::string>("folder")),
206  trackMinPt_(iConfig.getParameter<double>("trackMinimumPt")),
207  trackMaxBtlEta_(iConfig.getParameter<double>("trackMaximumBtlEta")),
208  trackMinEtlEta_(iConfig.getParameter<double>("trackMinimumEtlEta")),
209  trackMaxEtlEta_(iConfig.getParameter<double>("trackMaximumEtlEta")),
210  optionalPlots_(iConfig.getUntrackedParameter<bool>("optionalPlots")) {
211  GenRecTrackToken_ = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("inputTagG"));
212  RecTrackToken_ = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("inputTagT"));
213  RecVertexToken_ = consumes<std::vector<reco::Vertex>>(iConfig.getParameter<edm::InputTag>("inputTagV"));
214  HepMCProductToken_ = consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("inputTagH"));
216  consumes<TrackingParticleCollection>(iConfig.getParameter<edm::InputTag>("SimTag"));
218  consumes<reco::SimToRecoCollection>(iConfig.getParameter<edm::InputTag>("TPtoRecoTrackAssoc"));
220  consumes<reco::RecoToSimCollection>(iConfig.getParameter<edm::InputTag>("TPtoRecoTrackAssoc"));
221  btlSimHitsToken_ = consumes<CrossingFrame<PSimHit>>(iConfig.getParameter<edm::InputTag>("btlSimHits"));
222  etlSimHitsToken_ = consumes<CrossingFrame<PSimHit>>(iConfig.getParameter<edm::InputTag>("etlSimHits"));
223  trackAssocToken_ = consumes<edm::ValueMap<int>>(iConfig.getParameter<edm::InputTag>("trackAssocSrc"));
224  pathLengthToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("pathLengthSrc"));
225  tmtdToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("tmtd"));
226  SigmatmtdToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmatmtd"));
227  t0SrcToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("t0Src"));
228  Sigmat0SrcToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmat0Src"));
229  t0PidToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("t0PID"));
230  Sigmat0PidToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmat0PID"));
231  t0SafePidToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("t0SafePID"));
232  Sigmat0SafePidToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmat0SafePID"));
233  trackMVAQualToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("trackMVAQual"));
234  mtdgeoToken_ = esConsumes<MTDGeometry, MTDDigiGeometryRecord>();
235  mtdtopoToken_ = esConsumes<MTDTopology, MTDTopologyRcd>();
236  particleTableToken_ = esConsumes<HepPDT::ParticleDataTable, edm::DefaultRecord>();
237 }
238 
240 
241 // ------------ method called for each event ------------
243  using namespace edm;
244  using namespace geant_units::operators;
245  using namespace std;
246 
247  auto geometryHandle = iSetup.getTransientHandle(mtdgeoToken_);
248  const MTDGeometry* geom = geometryHandle.product();
249  auto topologyHandle = iSetup.getTransientHandle(mtdtopoToken_);
250  const MTDTopology* topology = topologyHandle.product();
251 
252  bool topo1Dis = false;
253  bool topo2Dis = false;
255  topo1Dis = true;
256  } else {
257  topo2Dis = true;
258  }
259 
260  auto GenRecTrackHandle = makeValid(iEvent.getHandle(GenRecTrackToken_));
261  auto RecVertexHandle = makeValid(iEvent.getHandle(RecVertexToken_));
262 
263  std::unordered_map<uint32_t, MTDHit> m_btlHits;
264  std::unordered_map<uint32_t, MTDHit> m_etlHits;
265  std::unordered_map<uint32_t, std::set<unsigned long int>> m_btlTrkPerCell;
266  std::unordered_map<uint32_t, std::set<unsigned long int>> m_etlTrkPerCell;
267 
268  const auto& tMtd = iEvent.get(tmtdToken_);
269  const auto& SigmatMtd = iEvent.get(SigmatmtdToken_);
270  const auto& t0Src = iEvent.get(t0SrcToken_);
271  const auto& Sigmat0Src = iEvent.get(Sigmat0SrcToken_);
272  const auto& t0Pid = iEvent.get(t0PidToken_);
273  const auto& Sigmat0Pid = iEvent.get(Sigmat0PidToken_);
274  const auto& t0Safe = iEvent.get(t0SafePidToken_);
275  const auto& Sigmat0Safe = iEvent.get(Sigmat0SafePidToken_);
276  const auto& mtdQualMVA = iEvent.get(trackMVAQualToken_);
277  const auto& trackAssoc = iEvent.get(trackAssocToken_);
278  const auto& pathLength = iEvent.get(pathLengthToken_);
279 
280  const auto& primRecoVtx = *(RecVertexHandle.product()->begin());
281 
282  // generator level information (HepMC format)
283  auto GenEventHandle = makeValid(iEvent.getHandle(HepMCProductToken_));
284  const HepMC::GenEvent* mc = GenEventHandle->GetEvent();
285  double zsim = convertMmToCm((*(mc->vertices_begin()))->position().z());
286  double tsim = (*(mc->vertices_begin()))->position().t() * CLHEP::mm / CLHEP::c_light;
287 
288  auto pdt = iSetup.getHandle(particleTableToken_);
289  const HepPDT::ParticleDataTable* pdTable = pdt.product();
290 
291  auto simToRecoH = makeValid(iEvent.getHandle(simToRecoAssociationToken_));
292  s2r_ = simToRecoH.product();
293 
294  auto recoToSimH = makeValid(iEvent.getHandle(recoToSimAssociationToken_));
295  r2s_ = recoToSimH.product();
296 
297  // find all signal event trackId corresponding to an MTD simHit
298  std::unordered_set<unsigned long int> mtdTrackId;
299 
300  std::unordered_set<unsigned long int> tpTrackId;
301  auto tpHandle = makeValid(iEvent.getHandle(trackingParticleCollectionToken_));
302  TrackingParticleCollection tpColl = *(tpHandle.product());
303  for (const auto& tp : tpColl) {
304  if (tp.eventId().bunchCrossing() == 0 && tp.eventId().event() == 0) {
305  if (!mvaTPSel(tp))
306  continue;
307  if (optionalPlots_) {
308  if (std::abs(tp.eta()) < trackMaxBtlEta_) {
310  } else if ((std::abs(tp.eta()) < trackMaxEtlEta_) && (std::abs(tp.eta()) > trackMinEtlEta_)) {
312  }
313  }
314  for (const auto& simTrk : tp.g4Tracks()) {
315  auto const thisTId = uniqueId(simTrk.trackId(), simTrk.eventId());
316  tpTrackId.insert(thisTId);
317  LogDebug("MtdTracksValidation") << "TP simTrack id : " << thisTId;
318  }
319  }
320  }
321 
322  //Fill maps with simhits accumulated per DetId
323 
324  auto btlSimHitsHandle = makeValid(iEvent.getHandle(btlSimHitsToken_));
325  MixCollection<PSimHit> btlSimHits(btlSimHitsHandle.product());
326  for (auto const& simHit : btlSimHits) {
327  if (simHit.tof() < 0 || simHit.tof() > 25.)
328  continue;
329  DetId id = simHit.detUnitId();
330  auto const thisHId = uniqueId(simHit.trackId(), simHit.eventId());
331  m_btlTrkPerCell[id.rawId()].insert(thisHId);
332  auto simHitIt = m_btlHits.emplace(id.rawId(), MTDHit()).first;
333  // --- Accumulate the energy (in MeV) of SIM hits in the same detector cell
334  (simHitIt->second).energy += convertUnitsTo(0.001_MeV, simHit.energyLoss());
335  }
336 
337  uint32_t hcount(0);
338  for (auto const& cell : m_btlTrkPerCell) {
339  bool foundAssocTP = false;
340  auto detId_key = cell.first;
341  for (auto const& simtrack : cell.second) {
342  if (tpTrackId.find(simtrack) != tpTrackId.end()) {
343  foundAssocTP = true;
344  mtdTrackId.insert(simtrack);
345  }
346  }
347  if (foundAssocTP == false) {
348  meUnassCrysEnergy_->Fill(log10(m_btlHits[detId_key].energy));
349  if (m_btlHits[detId_key].energy > depositBTLthreshold_) {
350  hcount++;
351  }
352  }
353  }
354  meUnassociatedDetId_->Fill(0.5, hcount);
355 
356  auto etlSimHitsHandle = makeValid(iEvent.getHandle(etlSimHitsToken_));
357  MixCollection<PSimHit> etlSimHits(etlSimHitsHandle.product());
358  for (auto const& simHit : etlSimHits) {
359  if (simHit.tof() < 0 || simHit.tof() > 25.) {
360  continue;
361  }
362  DetId id = simHit.detUnitId();
363  auto const thisHId = uniqueId(simHit.trackId(), simHit.eventId());
364  m_etlTrkPerCell[id.rawId()].insert(thisHId);
365  auto simHitIt = m_etlHits.emplace(id.rawId(), MTDHit()).first;
366  // --- Accumulate the energy (in MeV) of SIM hits in the same detector cell
367  (simHitIt->second).energy += convertUnitsTo(0.001_MeV, simHit.energyLoss());
368  }
369 
370  hcount = 0;
371  for (auto const& cell : m_etlTrkPerCell) {
372  bool foundAssocTP = false;
373  auto detId_key = cell.first;
374  for (auto const& simtrack : cell.second) {
375  if (tpTrackId.find(simtrack) != tpTrackId.end()) {
376  foundAssocTP = true;
377  mtdTrackId.insert(simtrack);
378  }
379  }
380  if (foundAssocTP == false) {
381  meUnassLgadsEnergy_->Fill(log10(m_etlHits[detId_key].energy));
382  if (m_etlHits[detId_key].energy > depositETLthreshold_) {
383  hcount++;
384  }
385  }
386  }
387  meUnassociatedDetId_->Fill(1.5, hcount);
388 
389  // Search for TP without associated hits and unassociated DetIds above threshold
390  if (optionalPlots_) {
391  for (const auto& tp : tpColl) {
392  if (tp.eventId().bunchCrossing() == 0 && tp.eventId().event() == 0) {
393  if (!mvaTPSel(tp)) {
394  continue;
395  }
396 
397  // test BTL crystals for association
398 
399  if (std::abs(tp.eta()) < trackMaxBtlEta_) {
400  bool tpIsAssoc = false;
401  bool goodCell = false;
402  for (auto const& cell : m_btlTrkPerCell) {
403  auto detId_key = cell.first;
404  if (m_btlHits[detId_key].energy < depositBTLthreshold_) {
405  continue;
406  }
407 
408  BTLDetId detId(detId_key);
410  const MTDGeomDet* thedet = geom->idToDet(geoId);
411  if (thedet == nullptr)
412  throw cms::Exception("MtdTracksValidation") << "GeographicalID: " << std::hex << geoId.rawId() << " ("
413  << detId.rawId() << ") is invalid!" << std::dec << std::endl;
414  const ProxyMTDTopology& topoproxy = static_cast<const ProxyMTDTopology&>(thedet->topology());
415  const RectangularMTDTopology& topo =
416  static_cast<const RectangularMTDTopology&>(topoproxy.specificTopology());
417 
418  Local3DPoint local_point(convertMmToCm(m_btlHits[detId_key].x),
419  convertMmToCm(m_btlHits[detId_key].y),
420  convertMmToCm(m_btlHits[detId_key].z));
421 
422  local_point =
423  topo.pixelToModuleLocalPoint(local_point, detId.row(topo.nrows()), detId.column(topo.nrows()));
424  const auto& global_point = thedet->toGlobal(local_point);
425 
426  if (std::abs(tp.eta() - global_point.eta()) > etaMatchCut_) {
427  continue;
428  }
429  goodCell = true;
430  for (auto const& simtrack : cell.second) {
431  for (auto const& TPsimtrack : tp.g4Tracks()) {
432  auto const testId = uniqueId(TPsimtrack.trackId(), TPsimtrack.eventId());
433  if (simtrack == testId) {
434  tpIsAssoc = true;
435  break;
436  }
437  }
438  }
439  } //cell Loop
440  if (!tpIsAssoc && goodCell) {
441  meUnassDeposit_->Fill(0.5);
442  }
443 
444  } else {
445  // test ETL LGADs for association
446  bool tpIsAssoc = false;
447  bool goodCell = false;
448  for (auto const& cell : m_etlTrkPerCell) {
449  auto detId_key = cell.first;
450  if (m_etlHits[detId_key].energy < depositETLthreshold_) {
451  continue;
452  }
453 
454  ETLDetId detId(detId_key);
455  DetId geoId = detId.geographicalId();
456  const MTDGeomDet* thedet = geom->idToDet(geoId);
457  if (thedet == nullptr)
458  throw cms::Exception("MtdTracksValidation") << "GeographicalID: " << std::hex << geoId.rawId() << " ("
459  << detId.rawId() << ") is invalid!" << std::dec << std::endl;
460 
461  Local3DPoint local_point(convertMmToCm(m_etlHits[detId_key].x),
462  convertMmToCm(m_etlHits[detId_key].y),
463  convertMmToCm(m_etlHits[detId_key].z));
464  const auto& global_point = thedet->toGlobal(local_point);
465 
466  if (std::abs(tp.eta() - global_point.eta()) > etaMatchCut_) {
467  continue;
468  }
469  goodCell = true;
470  for (auto const& simtrack : cell.second) {
471  for (auto const& TPsimtrack : tp.g4Tracks()) {
472  auto const testId = uniqueId(TPsimtrack.trackId(), TPsimtrack.eventId());
473  if (simtrack == testId) {
474  tpIsAssoc = true;
475  break;
476  }
477  }
478  }
479  } //cell Loop
480  if (!tpIsAssoc && goodCell) {
481  meUnassDeposit_->Fill(1.5);
482  }
483  } // tp BTL/ETL acceptance selection
484  }
485  } //tp Loop
486  } //optionalPlots
487 
488  unsigned int index = 0;
489 
490  // flag to select events with reco vertex close to true simulated primary vertex, or PV fake (particle guns)
491  const bool isGoodVtx = std::abs(primRecoVtx.z() - zsim) < deltaZcut_ || primRecoVtx.isFake();
492 
493  // --- Loop over all RECO tracks ---
494  for (const auto& trackGen : *GenRecTrackHandle) {
495  const reco::TrackRef trackref(iEvent.getHandle(GenRecTrackToken_), index);
496  index++;
497 
498  if (trackAssoc[trackref] == -1) {
499  LogInfo("mtdTracks") << "Extended track not associated";
500  continue;
501  }
502 
503  const reco::TrackRef mtdTrackref = reco::TrackRef(iEvent.getHandle(RecTrackToken_), trackAssoc[trackref]);
504  const reco::Track track = *mtdTrackref;
505 
506  bool isBTL = false;
507  bool twoETLdiscs = false;
508 
509  if (track.pt() >= trackMinPt_ && std::abs(track.eta()) <= trackMaxEtlEta_) {
510  meTracktmtd_->Fill(tMtd[trackref]);
511  if (std::round(SigmatMtd[trackref] - Sigmat0Pid[trackref]) != 0) {
512  LogWarning("mtdTracks")
513  << "TimeError associated to refitted track is different from TimeError stored in tofPID "
514  "sigmat0 ValueMap: this should not happen";
515  }
516 
517  meTrackt0Src_->Fill(t0Src[trackref]);
518  meTrackSigmat0Src_->Fill(Sigmat0Src[trackref]);
519 
520  meTrackt0Pid_->Fill(t0Pid[trackref]);
521  meTrackSigmat0Pid_->Fill(Sigmat0Pid[trackref]);
522  meTrackt0SafePid_->Fill(t0Safe[trackref]);
523  meTrackSigmat0SafePid_->Fill(Sigmat0Safe[trackref]);
524  meTrackMVAQual_->Fill(mtdQualMVA[trackref]);
525 
526  meTrackPathLenghtvsEta_->Fill(std::abs(track.eta()), pathLength[trackref]);
527 
528  if (std::abs(track.eta()) < trackMaxBtlEta_) {
529  // --- all BTL tracks (with and without hit in MTD) ---
533 
534  bool MTDBtl = false;
535  int numMTDBtlvalidhits = 0;
536  for (const auto hit : track.recHits()) {
537  if (hit->isValid() == false)
538  continue;
539  MTDDetId Hit = hit->geographicalId();
540  if ((Hit.det() == 6) && (Hit.subdetId() == 1) && (Hit.mtdSubDetector() == 1)) {
541  MTDBtl = true;
542  numMTDBtlvalidhits++;
543  }
544  }
545  meTrackNumHits_->Fill(numMTDBtlvalidhits);
546 
547  // --- keeping only tracks with last hit in MTD ---
548  if (MTDBtl == true) {
549  isBTL = true;
554  meBTLTrackPtRes_->Fill((trackGen.pt() - track.pt()) / trackGen.pt());
555  }
556  } //loop over (geometrical) BTL tracks
557 
558  else {
559  // --- all ETL tracks (with and without hit in MTD) ---
560  if ((track.eta() < -trackMinEtlEta_) && (track.eta() > -trackMaxEtlEta_)) {
561  meETLTrackEffEtaTot_[0]->Fill(track.eta());
562  meETLTrackEffPhiTot_[0]->Fill(track.phi());
563  meETLTrackEffPtTot_[0]->Fill(track.pt());
564  }
565 
566  if ((track.eta() > trackMinEtlEta_) && (track.eta() < trackMaxEtlEta_)) {
567  meETLTrackEffEtaTot_[1]->Fill(track.eta());
568  meETLTrackEffPhiTot_[1]->Fill(track.phi());
569  meETLTrackEffPtTot_[1]->Fill(track.pt());
570  }
571 
572  bool MTDEtlZnegD1 = false;
573  bool MTDEtlZnegD2 = false;
574  bool MTDEtlZposD1 = false;
575  bool MTDEtlZposD2 = false;
576  int numMTDEtlvalidhits = 0;
577  for (const auto hit : track.recHits()) {
578  if (hit->isValid() == false)
579  continue;
580  MTDDetId Hit = hit->geographicalId();
581  if ((Hit.det() == 6) && (Hit.subdetId() == 1) && (Hit.mtdSubDetector() == 2)) {
582  ETLDetId ETLHit = hit->geographicalId();
583 
584  if (topo2Dis) {
585  if ((ETLHit.zside() == -1) && (ETLHit.nDisc() == 1)) {
586  MTDEtlZnegD1 = true;
588  meETLTrackPtRes_->Fill((trackGen.pt() - track.pt()) / trackGen.pt());
589  numMTDEtlvalidhits++;
590  }
591  if ((ETLHit.zside() == -1) && (ETLHit.nDisc() == 2)) {
592  MTDEtlZnegD2 = true;
594  meETLTrackPtRes_->Fill((trackGen.pt() - track.pt()) / trackGen.pt());
595  numMTDEtlvalidhits++;
596  }
597  if ((ETLHit.zside() == 1) && (ETLHit.nDisc() == 1)) {
598  MTDEtlZposD1 = true;
600  meETLTrackPtRes_->Fill((trackGen.pt() - track.pt()) / trackGen.pt());
601  numMTDEtlvalidhits++;
602  }
603  if ((ETLHit.zside() == 1) && (ETLHit.nDisc() == 2)) {
604  MTDEtlZposD2 = true;
606  meETLTrackPtRes_->Fill((trackGen.pt() - track.pt()) / trackGen.pt());
607  numMTDEtlvalidhits++;
608  }
609  }
610 
611  if (topo1Dis) {
612  if (ETLHit.zside() == -1) {
613  MTDEtlZnegD1 = true;
615  numMTDEtlvalidhits++;
616  }
617  if (ETLHit.zside() == 1) {
618  MTDEtlZposD1 = true;
620  numMTDEtlvalidhits++;
621  }
622  }
623  }
624  }
625  meTrackNumHits_->Fill(-numMTDEtlvalidhits);
626 
627  // --- keeping only tracks with last hit in MTD ---
628  if ((track.eta() < -trackMinEtlEta_) && (track.eta() > -trackMaxEtlEta_)) {
629  twoETLdiscs = (MTDEtlZnegD1 == true) && (MTDEtlZnegD2 == true);
630  if ((MTDEtlZnegD1 == true) || (MTDEtlZnegD2 == true)) {
631  meETLTrackEffEtaMtd_[0]->Fill(track.eta());
632  meETLTrackEffPhiMtd_[0]->Fill(track.phi());
633  meETLTrackEffPtMtd_[0]->Fill(track.pt());
634  if (twoETLdiscs) {
635  meETLTrackEffEta2Mtd_[0]->Fill(track.eta());
636  meETLTrackEffPhi2Mtd_[0]->Fill(track.phi());
637  meETLTrackEffPt2Mtd_[0]->Fill(track.pt());
638  }
639  }
640  }
641  if ((track.eta() > trackMinEtlEta_) && (track.eta() < trackMaxEtlEta_)) {
642  twoETLdiscs = (MTDEtlZposD1 == true) && (MTDEtlZposD2 == true);
643  if ((MTDEtlZposD1 == true) || (MTDEtlZposD2 == true)) {
644  meETLTrackEffEtaMtd_[1]->Fill(track.eta());
645  meETLTrackEffPhiMtd_[1]->Fill(track.phi());
646  meETLTrackEffPtMtd_[1]->Fill(track.pt());
647  if (twoETLdiscs) {
648  meETLTrackEffEta2Mtd_[1]->Fill(track.eta());
649  meETLTrackEffPhi2Mtd_[1]->Fill(track.phi());
650  meETLTrackEffPt2Mtd_[1]->Fill(track.pt());
651  }
652  }
653  }
654  }
655  }
656 
657  if (isGoodVtx) {
658  bool noCrack = std::abs(trackGen.eta()) < trackMaxBtlEta_ || std::abs(trackGen.eta()) > trackMinEtlEta_;
659  const bool vtxFake = primRecoVtx.isFake();
660 
661  if (mvaRecSel(trackGen, primRecoVtx, t0Safe[trackref], Sigmat0Safe[trackref])) {
662  // reco-gen matching used for MVA quality flag
663 
664  if (noCrack) {
665  meMVATrackEffPtTot_->Fill(trackGen.pt());
666  }
667  meMVATrackEffEtaTot_->Fill(std::abs(trackGen.eta()));
668 
669  double dZ = trackGen.vz() - zsim;
670  double dT(-9999.);
671  double pullT(-9999.);
672  if (Sigmat0Safe[trackref] != -1.) {
673  dT = t0Safe[trackref] - tsim;
674  pullT = dT / Sigmat0Safe[trackref];
675  }
676  for (const auto& genP : mc->particle_range()) {
677  // select status 1 genParticles and match them to the reconstructed track
678 
679  float charge = pdTable->particle(HepPDT::ParticleID(genP->pdg_id())) != nullptr
680  ? pdTable->particle(HepPDT::ParticleID(genP->pdg_id()))->charge()
681  : 0.f;
682  if (mvaGenSel(*genP, charge)) {
683  if (mvaGenRecMatch(*genP, zsim, trackGen, vtxFake)) {
685  if (noCrack) {
686  meMVATrackMatchedEffPtTot_->Fill(trackGen.pt());
687  }
688  meMVATrackMatchedEffEtaTot_->Fill(std::abs(trackGen.eta()));
689  if (pullT > -9999.) {
690  meMVATrackResTot_->Fill(dT);
691  meMVATrackPullTot_->Fill(pullT);
692  if (noCrack) {
693  meMVATrackMatchedEffPtMtd_->Fill(trackGen.pt());
694  }
695  meMVATrackMatchedEffEtaMtd_->Fill(std::abs(trackGen.eta()));
696  }
697  break;
698  }
699  }
700  }
701 
702  // TrackingParticle based matching
703 
704  const reco::TrackBaseRef trkrefb(trackref);
705  auto tp_info = getMatchedTP(trkrefb, zsim);
706 
707  if (tp_info != nullptr && mvaTPSel(**tp_info)) {
708  const bool withMTD = tpWithMTD(**tp_info, mtdTrackId);
709  if (noCrack) {
710  meTrackMatchedTPEffPtTot_->Fill(trackGen.pt());
711  if (withMTD) {
712  meTrackMatchedTPmtdEffPtTot_->Fill(trackGen.pt());
713  }
714  }
715  meTrackMatchedTPEffEtaTot_->Fill(std::abs(trackGen.eta()));
716  if (withMTD) {
717  meTrackMatchedTPmtdEffEtaTot_->Fill(std::abs(trackGen.eta()));
718  }
719  if (pullT > -9999.) {
720  if (noCrack) {
721  meTrackMatchedTPEffPtMtd_->Fill(trackGen.pt());
722  if (isBTL || twoETLdiscs) {
723  meTrackMatchedTPEffPtEtl2Mtd_->Fill(trackGen.pt());
724  }
725  if (withMTD) {
726  meTrackMatchedTPmtdEffPtMtd_->Fill(trackGen.pt());
727  }
728  }
729  meTrackMatchedTPEffEtaMtd_->Fill(std::abs(trackGen.eta()));
730  if (isBTL || twoETLdiscs) {
731  meTrackMatchedTPEffEtaEtl2Mtd_->Fill(std::abs(trackGen.eta()));
732  }
733  if (withMTD) {
734  meTrackMatchedTPmtdEffEtaMtd_->Fill(std::abs(trackGen.eta()));
735  }
736  }
737  }
738  }
739  } // MC truth matich analysis for good PV
740  } //RECO tracks loop
741 }
742 
743 // ------------ method for histogram booking ------------
745  ibook.setCurrentFolder(folder_);
746 
747  // histogram booking
748  meBTLTrackRPTime_ = ibook.book1D("TrackBTLRPTime", "Track t0 with respect to R.P.;t0 [ns]", 100, -1, 3);
749  meBTLTrackEffEtaTot_ = ibook.book1D("TrackBTLEffEtaTot", "Track efficiency vs eta (Tot);#eta_{RECO}", 100, -1.6, 1.6);
751  ibook.book1D("TrackBTLEffPhiTot", "Track efficiency vs phi (Tot);#phi_{RECO} [rad]", 100, -3.2, 3.2);
752  meBTLTrackEffPtTot_ = ibook.book1D("TrackBTLEffPtTot", "Track efficiency vs pt (Tot);pt_{RECO} [GeV]", 50, 0, 10);
753  meBTLTrackEffEtaMtd_ = ibook.book1D("TrackBTLEffEtaMtd", "Track efficiency vs eta (Mtd);#eta_{RECO}", 100, -1.6, 1.6);
755  ibook.book1D("TrackBTLEffPhiMtd", "Track efficiency vs phi (Mtd);#phi_{RECO} [rad]", 100, -3.2, 3.2);
756  meBTLTrackEffPtMtd_ = ibook.book1D("TrackBTLEffPtMtd", "Track efficiency vs pt (Mtd);pt_{RECO} [GeV]", 50, 0, 10);
758  ibook.book1D("TrackBTLPtRes", "Track pT resolution ;pT_{Gentrack}-pT_{MTDtrack}/pT_{Gentrack} ", 100, -0.1, 0.1);
759  meETLTrackRPTime_ = ibook.book1D("TrackETLRPTime", "Track t0 with respect to R.P.;t0 [ns]", 100, -1, 3);
761  ibook.book1D("TrackETLEffEtaTotZneg", "Track efficiency vs eta (Tot) (-Z);#eta_{RECO}", 100, -3.2, -1.4);
763  ibook.book1D("TrackETLEffEtaTotZpos", "Track efficiency vs eta (Tot) (+Z);#eta_{RECO}", 100, 1.4, 3.2);
765  ibook.book1D("TrackETLEffPhiTotZneg", "Track efficiency vs phi (Tot) (-Z);#phi_{RECO} [rad]", 100, -3.2, 3.2);
767  ibook.book1D("TrackETLEffPhiTotZpos", "Track efficiency vs phi (Tot) (+Z);#phi_{RECO} [rad]", 100, -3.2, 3.2);
769  ibook.book1D("TrackETLEffPtTotZneg", "Track efficiency vs pt (Tot) (-Z);pt_{RECO} [GeV]", 50, 0, 10);
771  ibook.book1D("TrackETLEffPtTotZpos", "Track efficiency vs pt (Tot) (+Z);pt_{RECO} [GeV]", 50, 0, 10);
773  ibook.book1D("TrackETLEffEtaMtdZneg", "Track efficiency vs eta (Mtd) (-Z);#eta_{RECO}", 100, -3.2, -1.4);
775  ibook.book1D("TrackETLEffEtaMtdZpos", "Track efficiency vs eta (Mtd) (+Z);#eta_{RECO}", 100, 1.4, 3.2);
777  ibook.book1D("TrackETLEffPhiMtdZneg", "Track efficiency vs phi (Mtd) (-Z);#phi_{RECO} [rad]", 100, -3.2, 3.2);
779  ibook.book1D("TrackETLEffPhiMtdZpos", "Track efficiency vs phi (Mtd) (+Z);#phi_{RECO} [rad]", 100, -3.2, 3.2);
781  ibook.book1D("TrackETLEffPtMtdZneg", "Track efficiency vs pt (Mtd) (-Z);pt_{RECO} [GeV]", 50, 0, 10);
783  ibook.book1D("TrackETLEffPtMtdZpos", "Track efficiency vs pt (Mtd) (+Z);pt_{RECO} [GeV]", 50, 0, 10);
785  ibook.book1D("TrackETLEffEta2MtdZneg", "Track efficiency vs eta (Mtd 2 hit) (-Z);#eta_{RECO}", 100, -3.2, -1.4);
787  ibook.book1D("TrackETLEffEta2MtdZpos", "Track efficiency vs eta (Mtd 2 hit) (+Z);#eta_{RECO}", 100, 1.4, 3.2);
788  meETLTrackEffPhi2Mtd_[0] = ibook.book1D(
789  "TrackETLEffPhi2MtdZneg", "Track efficiency vs phi (Mtd 2 hit) (-Z);#phi_{RECO} [rad]", 100, -3.2, 3.2);
790  meETLTrackEffPhi2Mtd_[1] = ibook.book1D(
791  "TrackETLEffPhi2MtdZpos", "Track efficiency vs phi (Mtd 2 hit) (+Z);#phi_{RECO} [rad]", 100, -3.2, 3.2);
793  ibook.book1D("TrackETLEffPt2MtdZneg", "Track efficiency vs pt (Mtd 2 hit) (-Z);pt_{RECO} [GeV]", 50, 0, 10);
795  ibook.book1D("TrackETLEffPt2MtdZpos", "Track efficiency vs pt (Mtd 2 hit) (+Z);pt_{RECO} [GeV]", 50, 0, 10);
797  ibook.book1D("TrackETLPtRes", "Track pT resolution;pT_{Gentrack}-pT_{MTDtrack}/pT_{Gentrack} ", 100, -0.1, 0.1);
798 
799  meTracktmtd_ = ibook.book1D("Tracktmtd", "Track time from TrackExtenderWithMTD;tmtd [ns]", 150, 1, 16);
800  meTrackt0Src_ = ibook.book1D("Trackt0Src", "Track time from TrackExtenderWithMTD;t0Src [ns]", 100, -1.5, 1.5);
802  ibook.book1D("TrackSigmat0Src", "Time Error from TrackExtenderWithMTD; #sigma_{t0Src} [ns]", 100, 0, 0.1);
803 
804  meTrackt0Pid_ = ibook.book1D("Trackt0Pid", "Track t0 as stored in TofPid;t0 [ns]", 100, -1, 1);
805  meTrackSigmat0Pid_ = ibook.book1D("TrackSigmat0Pid", "Sigmat0 as stored in TofPid; #sigma_{t0} [ns]", 100, 0, 0.1);
806  meTrackt0SafePid_ = ibook.book1D("Trackt0SafePID", "Track t0 Safe as stored in TofPid;t0 [ns]", 100, -1, 1);
808  ibook.book1D("TrackSigmat0SafePID", "Sigmat0 Safe as stored in TofPid; #sigma_{t0} [ns]", 100, 0, 0.1);
809  meTrackNumHits_ = ibook.book1D("TrackNumHits", "Number of valid MTD hits per track ; Number of hits", 10, -5, 5);
810  meTrackMVAQual_ = ibook.book1D("TrackMVAQual", "Track MVA Quality as stored in Value Map ; MVAQual", 100, 0, 1);
812  "TrackPathLenghtvsEta", "MTD Track pathlength vs MTD track Eta;|#eta|;Pathlength", 100, 0, 3.2, 100.0, 400.0, "S");
813 
814  meMVATrackEffPtTot_ = ibook.book1D("MVAEffPtTot", "Pt of tracks associated to LV; track pt [GeV] ", 110, 0., 11.);
816  ibook.book1D("MVAMatchedEffPtTot", "Pt of tracks associated to LV matched to GEN; track pt [GeV] ", 110, 0., 11.);
818  "MVAMatchedEffPtMtd", "Pt of tracks associated to LV matched to GEN with time; track pt [GeV] ", 110, 0., 11.);
819 
821  ibook.book1D("MatchedTPEffPtTot", "Pt of tracks associated to LV matched to TP; track pt [GeV] ", 110, 0., 11.);
823  "MatchedTPEffPtMtd", "Pt of tracks associated to LV matched to TP with time; track pt [GeV] ", 110, 0., 11.);
825  ibook.book1D("MatchedTPEffPtEtl2Mtd",
826  "Pt of tracks associated to LV matched to TP with time, 2 ETL hits; track pt [GeV] ",
827  110,
828  0.,
829  11.);
830 
832  "MatchedTPmtdEffPtTot", "Pt of tracks associated to LV matched to TP-mtd hit; track pt [GeV] ", 110, 0., 11.);
834  ibook.book1D("MatchedTPmtdEffPtMtd",
835  "Pt of tracks associated to LV matched to TP-mtd hit with time; track pt [GeV] ",
836  110,
837  0.,
838  11.);
839 
840  meMVATrackEffEtaTot_ = ibook.book1D("MVAEffEtaTot", "Eta of tracks associated to LV; track eta ", 66, 0., 3.3);
842  ibook.book1D("MVAMatchedEffEtaTot", "Eta of tracks associated to LV matched to GEN; track eta ", 66, 0., 3.3);
844  "MVAMatchedEffEtaMtd", "Eta of tracks associated to LV matched to GEN with time; track eta ", 66, 0., 3.3);
845 
847  ibook.book1D("MatchedTPEffEtaTot", "Eta of tracks associated to LV matched to TP; track eta ", 66, 0., 3.3);
849  "MatchedTPEffEtaMtd", "Eta of tracks associated to LV matched to TP with time; track eta ", 66, 0., 3.3);
851  ibook.book1D("MatchedTPEffEtaEtl2Mtd",
852  "Eta of tracks associated to LV matched to TP with time, 2 ETL hits; track eta ",
853  66,
854  0.,
855  3.3);
856 
858  "MatchedTPmtdEffEtaTot", "Eta of tracks associated to LV matched to TP-mtd hit; track eta ", 66, 0., 3.3);
860  ibook.book1D("MatchedTPmtdEffEtaMtd",
861  "Eta of tracks associated to LV matched to TP-mtd hit with time; track eta ",
862  66,
863  0.,
864  3.3);
865 
866  meMVATrackResTot_ = ibook.book1D(
867  "MVATrackRes", "t_{rec} - t_{sim} for LV associated tracks; t_{rec} - t_{sim} [ns] ", 120, -0.15, 0.15);
869  ibook.book1D("MVATrackPull", "Pull for associated tracks; (t_{rec}-t_{sim})/#sigma_{t}", 50, -5., 5.);
871  "MVATrackZposResTot", "Z_{PCA} - Z_{sim} for associated tracks;Z_{PCA} - Z_{sim} [cm] ", 100, -0.1, 0.1);
872 
874  "UnassociatedDetId", "Number of MTD cell not associated to any TP per event", 2, 0., 2., 0., 100000., "S");
875  meNTrackingParticles_ = ibook.book1D("NTrackingParticles", "Total #Tracking particles", 2, 0, 2);
877  ibook.book1D("UnassDeposit",
878  "#Tracking particles with deposit over threshold in MTD cell, but with no cell associated to TP;",
879  2,
880  0,
881  2);
883  ibook.book1D("UnassCrysEnergy",
884  "Energy deposit in BTL crystal with no associated SimTracks;log_{10}(Energy [MeV]) ",
885  100,
886  -3.5,
887  1.5);
888  meUnassLgadsEnergy_ = ibook.book1D("UnassLgadsEnergy",
889  "Energy deposit in ETL LGADs with no associated SimTracks;log_{10}(Energy [MeV]) ",
890  100,
891  -3.5,
892  1.5);
893 }
894 
895 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
896 
899 
900  desc.add<std::string>("folder", "MTD/Tracks");
901  desc.add<edm::InputTag>("inputTagG", edm::InputTag("generalTracks"));
902  desc.add<edm::InputTag>("inputTagT", edm::InputTag("trackExtenderWithMTD"));
903  desc.add<edm::InputTag>("inputTagV", edm::InputTag("offlinePrimaryVertices4D"));
904  desc.add<edm::InputTag>("inputTagH", edm::InputTag("generatorSmeared"));
905  desc.add<edm::InputTag>("SimTag", edm::InputTag("mix", "MergedTrackTruth"));
906  desc.add<edm::InputTag>("TPtoRecoTrackAssoc", edm::InputTag("trackingParticleRecoTrackAsssociation"));
907  desc.add<edm::InputTag>("btlSimHits", edm::InputTag("mix", "g4SimHitsFastTimerHitsBarrel"));
908  desc.add<edm::InputTag>("etlSimHits", edm::InputTag("mix", "g4SimHitsFastTimerHitsEndcap"));
909  desc.add<edm::InputTag>("tmtd", edm::InputTag("trackExtenderWithMTD:generalTracktmtd"));
910  desc.add<edm::InputTag>("sigmatmtd", edm::InputTag("trackExtenderWithMTD:generalTracksigmatmtd"));
911  desc.add<edm::InputTag>("t0Src", edm::InputTag("trackExtenderWithMTD:generalTrackt0"));
912  desc.add<edm::InputTag>("sigmat0Src", edm::InputTag("trackExtenderWithMTD:generalTracksigmat0"));
913  desc.add<edm::InputTag>("trackAssocSrc", edm::InputTag("trackExtenderWithMTD:generalTrackassoc"))
914  ->setComment("Association between General and MTD Extended tracks");
915  desc.add<edm::InputTag>("pathLengthSrc", edm::InputTag("trackExtenderWithMTD:generalTrackPathLength"));
916  desc.add<edm::InputTag>("t0SafePID", edm::InputTag("tofPID:t0safe"));
917  desc.add<edm::InputTag>("sigmat0SafePID", edm::InputTag("tofPID:sigmat0safe"));
918  desc.add<edm::InputTag>("sigmat0PID", edm::InputTag("tofPID:sigmat0"));
919  desc.add<edm::InputTag>("t0PID", edm::InputTag("tofPID:t0"));
920  desc.add<edm::InputTag>("trackMVAQual", edm::InputTag("mtdTrackQualityMVA:mtdQualMVA"));
921  desc.add<double>("trackMinimumPt", 0.7); // [GeV]
922  desc.add<double>("trackMaximumBtlEta", 1.5);
923  desc.add<double>("trackMinimumEtlEta", 1.6);
924  desc.add<double>("trackMaximumEtlEta", 3.);
925  desc.addUntracked<bool>("optionalPlots", true);
926 
927  descriptions.add("mtdTracksValid", desc);
928 }
929 
930 const bool MtdTracksValidation::mvaGenSel(const HepMC::GenParticle& gp, const float& charge) {
931  bool match = false;
932  if (gp.status() != 1) {
933  return match;
934  }
935  match = charge != 0.f && gp.momentum().perp() > pTcut_ && std::abs(gp.momentum().eta()) < etacutGEN_;
936  return match;
937 }
938 
940  bool match = false;
941  if (tp.status() != 1) {
942  return match;
943  }
944  auto x_pv = tp.parentVertex()->position().x();
945  auto y_pv = tp.parentVertex()->position().y();
946  auto z_pv = tp.parentVertex()->position().z();
947 
948  auto r_pv = std::sqrt(x_pv * x_pv + y_pv * y_pv);
949 
950  match = tp.charge() != 0 && tp.pt() > pTcut_ && std::abs(tp.eta()) < etacutGEN_ && r_pv < rBTL_ && z_pv < zETL_;
951  return match;
952 }
953 
955  const reco::Vertex& vtx,
956  const double& t0,
957  const double& st0) {
958  bool match = false;
959  match = trk.pt() > pTcut_ && std::abs(trk.eta()) < etacutREC_ &&
960  (std::abs(trk.vz() - vtx.z()) <= deltaZcut_ || vtx.isFake());
961  if (st0 > 0.) {
962  match = match && std::abs(t0 - vtx.t()) < 3. * st0;
963  }
964  return match;
965 }
966 
968  const double& zsim,
969  const reco::TrackBase& trk,
970  const bool& vtxFake) {
971  bool match = false;
972  double dR = reco::deltaR(genP.momentum(), trk.momentum());
973  double genPT = genP.momentum().perp();
974  match = std::abs(genPT - trk.pt()) < trk.pt() * deltaPTcut_ && dR < deltaDRcut_ &&
975  (std::abs(trk.vz() - zsim) < deltaZcut_ || vtxFake);
976  return match;
977 }
978 
980  const double& zsim) {
981  auto found = r2s_->find(recoTrack);
982 
983  // reco track not matched to any TP
984  if (found == r2s_->end())
985  return nullptr;
986 
987  //matched TP equal to any TP associated to signal sim vertex
988  for (const auto& tp : found->val) {
989  if (tp.first->eventId().bunchCrossing() == 0 && tp.first->eventId().event() == 0 &&
990  std::abs(tp.first->parentVertex()->position().z() - zsim) < deltaZcut_) {
991  return &tp.first;
992  }
993  }
994 
995  // reco track not matched to any TP from vertex
996  return nullptr;
997 }
998 
1000  const std::unordered_set<unsigned long int>& trkList) {
1001  for (const auto& simTrk : tp.g4Tracks()) {
1002  for (const auto& mtdTrk : trkList) {
1003  if (uniqueId(simTrk.trackId(), simTrk.eventId()) == mtdTrk) {
1004  return true;
1005  }
1006  }
1007  }
1008  return false;
1009 }
1010 
MonitorElement * meETLTrackEffEtaTot_[2]
uint8_t geoId(const VFATFrame &frame)
retrieve the GEO information for this channel
edm::EDGetTokenT< edm::ValueMap< int > > trackAssocToken_
static constexpr double etacutGEN_
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_
MonitorElement * meETLTrackRPTime_
edm::EDGetTokenT< edm::ValueMap< float > > t0SafePidToken_
std::string folder_
edm::EDGetTokenT< edm::ValueMap< float > > t0SrcToken_
static constexpr double rBTL_
MonitorElement * meUnassociatedDetId_
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
MonitorElement * meMVATrackMatchedEffEtaMtd_
const bool mvaRecSel(const reco::TrackBase &, const reco::Vertex &, const double &, const double &)
virtual const Topology & topology() const
Definition: GeomDet.cc:67
virtual const PixelTopology & specificTopology() const
MonitorElement * meNTrackingParticles_
static constexpr double etacutREC_
static constexpr double deltaPTcut_
MonitorElement * meBTLTrackPtRes_
int row(unsigned nrows=kCrystalsPerModuleTdr) const
Definition: BTLDetId.h:73
edm::EDGetTokenT< std::vector< reco::Vertex > > RecVertexToken_
MonitorElement * meETLTrackEffPtTot_[2]
MonitorElement * meBTLTrackEffEtaMtd_
MonitorElement * meETLTrackEffPhi2Mtd_[2]
edm::EDGetTokenT< edm::ValueMap< float > > tmtdToken_
MonitorElement * meMVATrackMatchedEffPtTot_
MonitorElement * meMVATrackEffPtTot_
MonitorElement * meMVATrackEffEtaTot_
const reco::SimToRecoCollection * s2r_
ETLDetId geographicalId() const
Definition: ETLDetId.h:110
const unsigned long int uniqueId(const uint32_t x, const EncodedEventId &y)
LocalPoint pixelToModuleLocalPoint(const LocalPoint &plp, int row, int col) const
MonitorElement * meTrackMatchedTPmtdEffEtaTot_
const reco::RecoToSimCollection * r2s_
const_iterator find(const key_type &k) const
find element with specified reference key
const bool mvaGenSel(const HepMC::GenParticle &, const float &)
const edm::Ref< std::vector< TrackingParticle > > * getMatchedTP(const reco::TrackBaseRef &, const double &)
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 * meBTLTrackEffPhiTot_
int iEvent
Definition: GenABIO.cc:224
MonitorElement * meTrackMatchedTPmtdEffEtaMtd_
MonitorElement * meTrackMatchedTPEffEtaTot_
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_
MonitorElement * meUnassLgadsEnergy_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isETL(const double eta) const
#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_
MonitorElement * meTrackMVAQual_
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:151
MonitorElement * meETLTrackEffPhiTot_[2]
ETLDetId::EtlLayout etlLayoutFromTopoMode(const int &topoMode)
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
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)
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_
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:118
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
MonitorElement * meUnassDeposit_
HLT enums.
int column(unsigned nrows=kCrystalsPerModuleTdr) const
Definition: BTLDetId.h:78
int nDisc() const
Definition: ETLDetId.h:124
double a
Definition: hdecay.h:119
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:18
MonitorElement * meUnassCrysEnergy_
ESTransientHandle< T > getTransientHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:162
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
BTLDetId::CrysLayout crysLayoutFromTopoMode(const int &topoMode)
edm::ESGetToken< MTDGeometry, MTDDigiGeometryRecord > mtdgeoToken_
const bool tpWithMTD(const TrackingParticle &, const std::unordered_set< unsigned long int > &)
std::vector< TrackingParticle > TrackingParticleCollection
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_
BTLDetId geographicalId(CrysLayout lay) const
Definition: BTLDetId.cc:3
Definition: Run.h:45
edm::EDGetTokenT< edm::ValueMap< float > > Sigmat0SafePidToken_
MtdTracksValidation(const edm::ParameterSet &)
#define LogDebug(id)
edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > particleTableToken_