Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 10:57:26

0001 #include "FWCore/ServiceRegistry/interface/Service.h"
0002 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0003 #include "DataFormats/MuonReco/interface/Muon.h"
0004 #include "DataFormats/MuonReco/interface/MuonFwd.h"
0005 #include "DataFormats/MuonReco/interface/MuonSelectors.h"
0006 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0007 #include "DataFormats/Math/interface/deltaR.h"
0008 #include "DataFormats/Math/interface/deltaPhi.h"
0009 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
0010 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0011 #include "DataFormats/EgammaCandidates/interface/Electron.h"
0012 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
0013 #include "DataFormats/Common/interface/Handle.h"
0014 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0015 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0016 #include "DataFormats/TrackReco/interface/Track.h"
0017 #include "DataFormats/VertexReco/interface/Vertex.h"
0018 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0019 #include "DQM/TrackingMonitorSource/interface/TrackTypeMonitor.h"
0020 
0021 #include "TFile.h"
0022 #include "TH1.h"
0023 #include "TMath.h"
0024 #include "TPRegexp.h"
0025 
0026 template <class T>
0027 class PtComparator {
0028 public:
0029   bool operator()(const T* a, const T* b) const { return (a->pt() > b->pt()); }
0030 };
0031 // -----------------------------
0032 // constructors and destructor
0033 // -----------------------------
0034 TrackTypeMonitor::TrackTypeMonitor(const edm::ParameterSet& ps)
0035     : parameters_(ps),
0036       moduleName_(parameters_.getUntrackedParameter<std::string>("ModuleName", "TrackTypeMonitor")),
0037       folderName_(parameters_.getUntrackedParameter<std::string>("FolderName", "highPurityTracks")),
0038       verbose_(parameters_.getUntrackedParameter<bool>("verbose", false)),
0039       muonTag_(ps.getUntrackedParameter<edm::InputTag>("muonInputTag", edm::InputTag("muons"))),
0040       electronTag_(ps.getUntrackedParameter<edm::InputTag>("electronInputTag", edm::InputTag("gedGsfElectrons"))),
0041       trackTag_(parameters_.getUntrackedParameter<edm::InputTag>("trackInputTag", edm::InputTag("generalTracks"))),
0042       bsTag_(parameters_.getUntrackedParameter<edm::InputTag>("offlineBeamSpot", edm::InputTag("offlineBeamSpot"))),
0043       vertexTag_(
0044           parameters_.getUntrackedParameter<edm::InputTag>("vertexTag", edm::InputTag("offlinePrimaryVertices"))),
0045       electronToken_(consumes<reco::GsfElectronCollection>(electronTag_)),
0046       muonToken_(consumes<reco::MuonCollection>(muonTag_)),
0047       trackToken_(consumes<reco::TrackCollection>(trackTag_)),
0048       bsToken_(consumes<reco::BeamSpot>(bsTag_)),
0049       vertexToken_(consumes<reco::VertexCollection>(vertexTag_)),
0050       trackQuality_(parameters_.getUntrackedParameter<std::string>("trackQuality", "highPurity")) {
0051   trackEtaHList_.clear();
0052   trackPhiHList_.clear();
0053   trackPHList_.clear();
0054   trackPtHList_.clear();
0055   trackPterrHList_.clear();
0056   trackqOverpHList_.clear();
0057   trackChi2bynDOFHList_.clear();
0058   nTracksHList_.clear();
0059 }
0060 void TrackTypeMonitor::bookHistograms(DQMStore::IBooker& iBook, edm::Run const& iRun, edm::EventSetup const& iSetup) {
0061   edm::ParameterSet TrackEtaHistoPar = parameters_.getParameter<edm::ParameterSet>("TrackEtaPar");
0062   edm::ParameterSet TrackPhiHistoPar = parameters_.getParameter<edm::ParameterSet>("TrackPhiPar");
0063   edm::ParameterSet TrackPHistoPar = parameters_.getParameter<edm::ParameterSet>("TrackPPar");
0064   edm::ParameterSet TrackPtHistoPar = parameters_.getParameter<edm::ParameterSet>("TrackPtPar");
0065   edm::ParameterSet TrackPterrHistoPar = parameters_.getParameter<edm::ParameterSet>("TrackPterrPar");
0066   edm::ParameterSet TrackqOverpHistoPar = parameters_.getParameter<edm::ParameterSet>("TrackqOverpPar");
0067   edm::ParameterSet TrackdzHistoPar = parameters_.getParameter<edm::ParameterSet>("TrackdzPar");
0068   edm::ParameterSet TrackChi2bynDOFHistoPar = parameters_.getParameter<edm::ParameterSet>("TrackChi2bynDOFPar");
0069   edm::ParameterSet nTracksHistoPar = parameters_.getParameter<edm::ParameterSet>("nTracksPar");
0070 
0071   std::string currentFolder = moduleName_ + "/" + folderName_;
0072   iBook.setCurrentFolder(currentFolder);
0073 
0074   trackEtaHList_.push_back(iBook.book1D("trackEtaIso",
0075                                         "Isolated Track Eta",
0076                                         TrackEtaHistoPar.getParameter<int32_t>("Xbins"),
0077                                         TrackEtaHistoPar.getParameter<double>("Xmin"),
0078                                         TrackEtaHistoPar.getParameter<double>("Xmax")));
0079   trackEtaHList_.push_back(iBook.book1D("trackEtaNoIso",
0080                                         "NonIsolated Track Eta",
0081                                         TrackEtaHistoPar.getParameter<int32_t>("Xbins"),
0082                                         TrackEtaHistoPar.getParameter<double>("Xmin"),
0083                                         TrackEtaHistoPar.getParameter<double>("Xmax")));
0084   trackEtaHList_.push_back(iBook.book1D("trackEtaUL",
0085                                         "Underlying Track Eta",
0086                                         TrackEtaHistoPar.getParameter<int32_t>("Xbins"),
0087                                         TrackEtaHistoPar.getParameter<double>("Xmin"),
0088                                         TrackEtaHistoPar.getParameter<double>("Xmax")));
0089   trackEtaHList_.push_back(iBook.book1D("trackEtaALL",
0090                                         "All Track Eta",
0091                                         TrackEtaHistoPar.getParameter<int32_t>("Xbins"),
0092                                         TrackEtaHistoPar.getParameter<double>("Xmin"),
0093                                         TrackEtaHistoPar.getParameter<double>("Xmax")));
0094 
0095   trackPhiHList_.push_back(iBook.book1D("trackPhiIso",
0096                                         "Isolated Track Phi",
0097                                         TrackPhiHistoPar.getParameter<int32_t>("Xbins"),
0098                                         TrackPhiHistoPar.getParameter<double>("Xmin"),
0099                                         TrackPhiHistoPar.getParameter<double>("Xmax")));
0100   trackPhiHList_.push_back(iBook.book1D("trackPhiNonIso",
0101                                         "NonIsolated Track Phi",
0102                                         TrackPhiHistoPar.getParameter<int32_t>("Xbins"),
0103                                         TrackPhiHistoPar.getParameter<double>("Xmin"),
0104                                         TrackPhiHistoPar.getParameter<double>("Xmax")));
0105   trackPhiHList_.push_back(iBook.book1D("trackPhiUL",
0106                                         "Underlying Track Phi",
0107                                         TrackPhiHistoPar.getParameter<int32_t>("Xbins"),
0108                                         TrackPhiHistoPar.getParameter<double>("Xmin"),
0109                                         TrackPhiHistoPar.getParameter<double>("Xmax")));
0110   trackPhiHList_.push_back(iBook.book1D("trackPhiALL",
0111                                         "All Track Phi",
0112                                         TrackPhiHistoPar.getParameter<int32_t>("Xbins"),
0113                                         TrackPhiHistoPar.getParameter<double>("Xmin"),
0114                                         TrackPhiHistoPar.getParameter<double>("Xmax")));
0115 
0116   trackPHList_.push_back(iBook.book1D("trackPIso",
0117                                       "Isolated Track P",
0118                                       TrackPHistoPar.getParameter<int32_t>("Xbins"),
0119                                       TrackPHistoPar.getParameter<double>("Xmin"),
0120                                       TrackPHistoPar.getParameter<double>("Xmax")));
0121   trackPHList_.push_back(iBook.book1D("trackPNonIso",
0122                                       "NonIsolated Track P",
0123                                       TrackPHistoPar.getParameter<int32_t>("Xbins"),
0124                                       TrackPHistoPar.getParameter<double>("Xmin"),
0125                                       TrackPHistoPar.getParameter<double>("Xmax")));
0126   trackPHList_.push_back(iBook.book1D("trackPUL",
0127                                       "Underlying Track P",
0128                                       TrackPHistoPar.getParameter<int32_t>("Xbins"),
0129                                       TrackPHistoPar.getParameter<double>("Xmin"),
0130                                       TrackPHistoPar.getParameter<double>("Xmax")));
0131   trackPHList_.push_back(iBook.book1D("trackPALL",
0132                                       "All Track P",
0133                                       TrackPHistoPar.getParameter<int32_t>("Xbins"),
0134                                       TrackPHistoPar.getParameter<double>("Xmin"),
0135                                       TrackPHistoPar.getParameter<double>("Xmax")));
0136 
0137   trackPtHList_.push_back(iBook.book1D("trackPtIsolated",
0138                                        "Isolated Track Pt",
0139                                        TrackPtHistoPar.getParameter<int32_t>("Xbins"),
0140                                        TrackPtHistoPar.getParameter<double>("Xmin"),
0141                                        TrackPtHistoPar.getParameter<double>("Xmax")));
0142   trackPtHList_.push_back(iBook.book1D("trackPtNonIsolated",
0143                                        "NonIsolated Track Pt",
0144                                        TrackPtHistoPar.getParameter<int32_t>("Xbins"),
0145                                        TrackPtHistoPar.getParameter<double>("Xmin"),
0146                                        TrackPtHistoPar.getParameter<double>("Xmax")));
0147   trackPtHList_.push_back(iBook.book1D("trackPtUL",
0148                                        "Underlying Track Pt",
0149                                        TrackPtHistoPar.getParameter<int32_t>("Xbins"),
0150                                        TrackPtHistoPar.getParameter<double>("Xmin"),
0151                                        TrackPtHistoPar.getParameter<double>("Xmax")));
0152   trackPtHList_.push_back(iBook.book1D("trackPtALL",
0153                                        "All Track Pt",
0154                                        TrackPtHistoPar.getParameter<int32_t>("Xbins"),
0155                                        TrackPtHistoPar.getParameter<double>("Xmin"),
0156                                        TrackPtHistoPar.getParameter<double>("Xmax")));
0157 
0158   trackPterrHList_.push_back(iBook.book1D("trackPterrIsolated",
0159                                           "Isolated Track Pterr",
0160                                           TrackPterrHistoPar.getParameter<int32_t>("Xbins"),
0161                                           TrackPterrHistoPar.getParameter<double>("Xmin"),
0162                                           TrackPterrHistoPar.getParameter<double>("Xmax")));
0163   trackPterrHList_.push_back(iBook.book1D("trackPterrNonIsolated",
0164                                           "NonIsolated Track Pterr",
0165                                           TrackPterrHistoPar.getParameter<int32_t>("Xbins"),
0166                                           TrackPterrHistoPar.getParameter<double>("Xmin"),
0167                                           TrackPterrHistoPar.getParameter<double>("Xmax")));
0168   trackPterrHList_.push_back(iBook.book1D("trackPterrUL",
0169                                           "Underlying Track Pterr",
0170                                           TrackPterrHistoPar.getParameter<int32_t>("Xbins"),
0171                                           TrackPterrHistoPar.getParameter<double>("Xmin"),
0172                                           TrackPterrHistoPar.getParameter<double>("Xmax")));
0173   trackPterrHList_.push_back(iBook.book1D("trackPterrALL",
0174                                           "All Track Pterr",
0175                                           TrackPterrHistoPar.getParameter<int32_t>("Xbins"),
0176                                           TrackPterrHistoPar.getParameter<double>("Xmin"),
0177                                           TrackPterrHistoPar.getParameter<double>("Xmax")));
0178 
0179   trackqOverpHList_.push_back(iBook.book1D("trackqOverpIsolated",
0180                                            "Isolated Track qOverp",
0181                                            TrackqOverpHistoPar.getParameter<int32_t>("Xbins"),
0182                                            TrackqOverpHistoPar.getParameter<double>("Xmin"),
0183                                            TrackqOverpHistoPar.getParameter<double>("Xmax")));
0184   trackqOverpHList_.push_back(iBook.book1D("trackqOverpNonIsolated",
0185                                            "NonIsolated Track qOverp",
0186                                            TrackqOverpHistoPar.getParameter<int32_t>("Xbins"),
0187                                            TrackqOverpHistoPar.getParameter<double>("Xmin"),
0188                                            TrackqOverpHistoPar.getParameter<double>("Xmax")));
0189   trackqOverpHList_.push_back(iBook.book1D("trackqOverpUL",
0190                                            "Underlying Track qOverp",
0191                                            TrackqOverpHistoPar.getParameter<int32_t>("Xbins"),
0192                                            TrackqOverpHistoPar.getParameter<double>("Xmin"),
0193                                            TrackqOverpHistoPar.getParameter<double>("Xmax")));
0194   trackqOverpHList_.push_back(iBook.book1D("trackqOverpALL",
0195                                            "All Track qOverp",
0196                                            TrackqOverpHistoPar.getParameter<int32_t>("Xbins"),
0197                                            TrackqOverpHistoPar.getParameter<double>("Xmin"),
0198                                            TrackqOverpHistoPar.getParameter<double>("Xmax")));
0199 
0200   trackdzHList_.push_back(iBook.book1D("trackdzIsolated",
0201                                        "Isolated Track dz",
0202                                        TrackdzHistoPar.getParameter<int32_t>("Xbins"),
0203                                        TrackdzHistoPar.getParameter<double>("Xmin"),
0204                                        TrackdzHistoPar.getParameter<double>("Xmax")));
0205   trackdzHList_.push_back(iBook.book1D("trackdzNonIsolated",
0206                                        "NonIsolated Track dz",
0207                                        TrackdzHistoPar.getParameter<int32_t>("Xbins"),
0208                                        TrackdzHistoPar.getParameter<double>("Xmin"),
0209                                        TrackdzHistoPar.getParameter<double>("Xmax")));
0210   trackdzHList_.push_back(iBook.book1D("trackdzUL",
0211                                        "Underlying Track dz",
0212                                        TrackdzHistoPar.getParameter<int32_t>("Xbins"),
0213                                        TrackdzHistoPar.getParameter<double>("Xmin"),
0214                                        TrackdzHistoPar.getParameter<double>("Xmax")));
0215   trackdzHList_.push_back(iBook.book1D("trackdzALL",
0216                                        "All Track dz",
0217                                        TrackdzHistoPar.getParameter<int32_t>("Xbins"),
0218                                        TrackdzHistoPar.getParameter<double>("Xmin"),
0219                                        TrackdzHistoPar.getParameter<double>("Xmax")));
0220 
0221   trackChi2bynDOFHList_.push_back(iBook.book1D("trackChi2bynDOFIsolated",
0222                                                "Isolated Track Chi2bynDOF",
0223                                                TrackChi2bynDOFHistoPar.getParameter<int32_t>("Xbins"),
0224                                                TrackChi2bynDOFHistoPar.getParameter<double>("Xmin"),
0225                                                TrackChi2bynDOFHistoPar.getParameter<double>("Xmax")));
0226   trackChi2bynDOFHList_.push_back(iBook.book1D("trackChi2bynDOFNonIsolated",
0227                                                "NonIsolated Track Chi2bynDOF",
0228                                                TrackChi2bynDOFHistoPar.getParameter<int32_t>("Xbins"),
0229                                                TrackChi2bynDOFHistoPar.getParameter<double>("Xmin"),
0230                                                TrackChi2bynDOFHistoPar.getParameter<double>("Xmax")));
0231   trackChi2bynDOFHList_.push_back(iBook.book1D("trackChi2bynDOFUL",
0232                                                "Underlying Track Chi2bynDOF",
0233                                                TrackChi2bynDOFHistoPar.getParameter<int32_t>("Xbins"),
0234                                                TrackChi2bynDOFHistoPar.getParameter<double>("Xmin"),
0235                                                TrackChi2bynDOFHistoPar.getParameter<double>("Xmax")));
0236   trackChi2bynDOFHList_.push_back(iBook.book1D("trackChi2bynDOFAll",
0237                                                "All Track Chi2bynDOF",
0238                                                TrackChi2bynDOFHistoPar.getParameter<int32_t>("Xbins"),
0239                                                TrackChi2bynDOFHistoPar.getParameter<double>("Xmin"),
0240                                                TrackChi2bynDOFHistoPar.getParameter<double>("Xmax")));
0241 
0242   nTracksHList_.push_back(iBook.book1D("nTracksIsolated",
0243                                        "Isolated Track nTracks",
0244                                        nTracksHistoPar.getParameter<int32_t>("Xbins"),
0245                                        nTracksHistoPar.getParameter<double>("Xmin"),
0246                                        nTracksHistoPar.getParameter<double>("Xmax")));
0247   nTracksHList_.push_back(iBook.book1D("nTracksNonIsolated",
0248                                        "NonIsolated Track nTracks",
0249                                        nTracksHistoPar.getParameter<int32_t>("Xbins"),
0250                                        nTracksHistoPar.getParameter<double>("Xmin"),
0251                                        nTracksHistoPar.getParameter<double>("Xmax")));
0252   nTracksHList_.push_back(iBook.book1D("nTracksUL",
0253                                        "Underlying Track nTracks",
0254                                        nTracksHistoPar.getParameter<int32_t>("Xbins"),
0255                                        nTracksHistoPar.getParameter<double>("Xmin"),
0256                                        nTracksHistoPar.getParameter<double>("Xmax")));
0257   nTracksHList_.push_back(iBook.book1D("nTracksAll",
0258                                        "All Track nTracks",
0259                                        nTracksHistoPar.getParameter<int32_t>("Xbins"),
0260                                        nTracksHistoPar.getParameter<double>("Xmin"),
0261                                        nTracksHistoPar.getParameter<double>("Xmax")));
0262 
0263   hcounterH_ = iBook.book1D("hcounter", "hcounter", 7, -0.5, 6.5);
0264   dphiH_ = iBook.book1D("dphi", "dphi", 100, 0, 7);
0265   drH_ = iBook.book1D("dr", "dr", 100, 0, 6);
0266 }
0267 void TrackTypeMonitor::analyze(edm::Event const& iEvent, edm::EventSetup const& iSetup) {
0268   // read the beam spot
0269   edm::Handle<reco::BeamSpot> beamSpot;
0270   iEvent.getByToken(bsToken_, beamSpot);
0271 
0272   std::vector<const reco::Track*> isoTrkList;
0273 
0274   // muons
0275   edm::Handle<reco::MuonCollection> muonColl;
0276   iEvent.getByToken(muonToken_, muonColl);
0277 
0278   if (muonColl.isValid()) {
0279     for (auto const& muo : *muonColl) {
0280       if (muo.isGlobalMuon() && muo.isPFMuon() && std::fabs(muo.eta()) <= 2.1 && muo.pt() > 1) {
0281         reco::TrackRef gtrkref = muo.globalTrack();
0282         if (!gtrkref.isNonnull())
0283           continue;
0284         const reco::Track* gtk = &(*gtrkref);
0285         double chi2 = gtk->chi2();
0286         double ndof = gtk->ndof();
0287         double chbyndof = (ndof > 0) ? chi2 / ndof : 0;
0288 
0289         const reco::HitPattern& hitp = gtk->hitPattern();
0290         int nPixelHits = hitp.numberOfValidPixelHits();
0291         int nStripHits = hitp.numberOfValidStripHits();
0292 
0293         reco::TrackRef itrkref = muo.innerTrack();  // tracker segment only
0294         if (!itrkref.isNonnull())
0295           continue;
0296         const reco::Track* tk = &(*itrkref);
0297         double trkd0 = tk->d0();
0298         double trkdz = tk->dz();
0299         if (beamSpot.isValid()) {
0300           trkd0 = -(tk->dxy(beamSpot->position()));
0301           trkdz = tk->dz(beamSpot->position());
0302         }
0303 
0304         // Hits/section in the muon chamber
0305         int nChambers = muo.numberOfChambers();
0306         int nMatches = muo.numberOfMatches();
0307         int nMatchedStations = muo.numberOfMatchedStations();
0308 
0309         // PF Isolation
0310         const reco::MuonPFIsolation& pfIso04 = muo.pfIsolationR04();
0311         double absiso = pfIso04.sumChargedParticlePt +
0312                         std::max(0.0, pfIso04.sumNeutralHadronEt + pfIso04.sumPhotonEt - 0.5 * pfIso04.sumPUPt);
0313         if (chbyndof < 10 && std::fabs(trkd0) < 0.02 && std::fabs(trkdz) < 20 && nPixelHits > 1 && nStripHits > 8 &&
0314             nChambers > 2 && nMatches > 2 && nMatchedStations > 2 && absiso / muo.pt() < 0.3) {
0315           isoTrkList.push_back(gtk);
0316           fillHistograms(*gtk, 0);
0317         }
0318       }
0319     }
0320   }
0321 
0322   // electrons
0323   edm::Handle<reco::GsfElectronCollection> electronColl;
0324   iEvent.getByToken(electronToken_, electronColl);
0325 
0326   if (electronColl.isValid()) {
0327     for (auto const& ele : *electronColl) {
0328       if (!ele.ecalDriven())
0329         continue;
0330       if (ele.pt() < 5)
0331         continue;
0332       hcounterH_->Fill(0);
0333 
0334       double hOverE = ele.hadronicOverEm();
0335       double sigmaee = ele.sigmaIetaIeta();
0336       double deltaPhiIn = ele.deltaPhiSuperClusterTrackAtVtx();
0337       double deltaEtaIn = ele.deltaEtaSuperClusterTrackAtVtx();
0338 
0339       if (ele.isEB()) {
0340         if (std::fabs(deltaPhiIn) >= .15 && std::fabs(deltaEtaIn) >= .007 && hOverE >= .12 && sigmaee >= .01)
0341           continue;
0342       } else if (ele.isEE()) {
0343         if (std::fabs(deltaPhiIn) >= .10 && std::fabs(deltaEtaIn) >= .009 && hOverE >= .10 && sigmaee >= .03)
0344           continue;
0345       }
0346       hcounterH_->Fill(1);
0347 
0348       reco::GsfTrackRef gsftrk = ele.gsfTrack();
0349       if (!gsftrk.isNonnull())
0350         continue;
0351       const reco::GsfTrack* trk = &(*gsftrk);
0352       double trkd0 = trk->d0();
0353       double trkdz = trk->dz();
0354       if (beamSpot.isValid()) {
0355         trkd0 = -(trk->dxy(beamSpot->position()));
0356         trkdz = trk->dz(beamSpot->position());
0357       }
0358       double chi2 = trk->chi2();
0359       double ndof = trk->ndof();
0360       double chbyndof = (ndof > 0) ? chi2 / ndof : 0;
0361       if (chbyndof >= 10)
0362         continue;
0363       hcounterH_->Fill(2);
0364 
0365       if (std::fabs(trkd0) >= 0.02 || std::fabs(trkdz) >= 20)
0366         continue;
0367       hcounterH_->Fill(3);
0368 
0369       const reco::HitPattern& hitp = trk->hitPattern();
0370       int nPixelHits = hitp.numberOfValidPixelHits();
0371       if (nPixelHits < 1)
0372         continue;
0373       hcounterH_->Fill(4);
0374 
0375       int nStripHits = hitp.numberOfValidStripHits();
0376       if (nStripHits < 8)
0377         continue;
0378       hcounterH_->Fill(5);
0379 
0380       reco::GsfElectron::PflowIsolationVariables pfIso = ele.pfIsolationVariables();
0381       float absiso =
0382           pfIso.sumChargedHadronPt + std::max(0.0, pfIso.sumNeutralHadronEt + pfIso.sumPhotonEt - 0.5 * pfIso.sumPUPt);
0383       float eiso = absiso / ele.pt();
0384       if (eiso > 0.2)
0385         continue;
0386       hcounterH_->Fill(6);
0387 
0388       isoTrkList.push_back(trk);
0389       fillHistograms(*trk, 0);
0390     }
0391   }
0392   nTracksHList_.at(0)->Fill(isoTrkList.size());
0393 
0394   // Read track collection
0395   edm::Handle<reco::TrackCollection> tracks;
0396   iEvent.getByToken(trackToken_, tracks);
0397 
0398   // Read vertex collection
0399   edm::Handle<reco::VertexCollection> vertexColl;
0400   iEvent.getByToken(vertexToken_, vertexColl);
0401 
0402   if (tracks.isValid()) {
0403     std::vector<const reco::Track*> nisoTrkList;
0404     const reco::Vertex& vit = vertexColl->front();  // Highest sumPt vertex
0405     reco::Track::TrackQuality quality = reco::Track::qualityByName(trackQuality_);
0406 
0407     int nmatch = 0;
0408     int nTracks = 0, nMBTracks = 0;
0409     for (auto const& track : *tracks) {
0410       if (!track.quality(quality))
0411         continue;
0412       ++nTracks;
0413       fillHistograms(track, 3);
0414 
0415       // now classify primary and underlying tracks using vertex information
0416       double dxy = track.dxy(vit.position());
0417       double dz = track.dz(vit.position());
0418       if (std::fabs(dxy) < 0.02 && std::fabs(dz) < 20) {  // primary tracks
0419                                                           // remove tracks associated to the identified particles
0420         double drmin = 999;
0421         for (auto const& tk : isoTrkList) {
0422           double dr = deltaR(track.eta(), track.phi(), tk->eta(), tk->phi());
0423           if (dr < drmin)
0424             drmin = dr;
0425         }
0426         if (drmin < 0.01) {
0427           if (verbose_) {
0428             edm::LogInfo("TrackTypeMonitor") << " Match: " << ++nmatch << " drmin = " << drmin;
0429             edm::LogInfo("TrackTypeMonitor")
0430                 << " Track Pt: " << track.pt() << " Eta: " << track.eta() << " Phi: " << track.phi();
0431             for (auto const& isotk : isoTrkList) {
0432               edm::LogInfo("TrackTypeMonitor")
0433                   << " Lepton Pt: " << isotk->pt() << " Eta: " << isotk->eta() << " Phi: " << isotk->phi();
0434             }
0435           }
0436           continue;
0437         }
0438 
0439         fillHistograms(track, 1);
0440         nisoTrkList.push_back(&track);
0441       } else {
0442         ++nMBTracks;
0443         fillHistograms(track, 2);  //non-primary tracks
0444       }
0445     }
0446     nTracksHList_.at(1)->Fill(nisoTrkList.size());
0447     nTracksHList_.at(2)->Fill(nMBTracks);
0448     nTracksHList_.at(3)->Fill(nTracks);
0449 
0450     std::sort(nisoTrkList.begin(), nisoTrkList.end(), PtComparator<reco::Track>());
0451     const reco::Track* isoTrk = isoTrkList.at(0);
0452     for (auto const& obj : nisoTrkList) {
0453       if (obj->pt() > 5) {
0454         double dphi = deltaPhi(isoTrk->phi(), obj->phi());
0455         dphiH_->Fill(std::fabs(dphi));
0456 
0457         double dr = deltaR(isoTrk->eta(), isoTrk->phi(), obj->eta(), obj->phi());
0458         drH_->Fill(dr);
0459       }
0460     }
0461   }
0462 }
0463 void TrackTypeMonitor::fillHistograms(const reco::Track& track, int indx) {
0464   if (indx >= 0 && indx < static_cast<int>(trackEtaHList_.size()))
0465     trackEtaHList_.at(indx)->Fill(track.eta());
0466 
0467   if (indx >= 0 && indx < static_cast<int>(trackPhiHList_.size()))
0468     trackPhiHList_.at(indx)->Fill(track.phi());
0469 
0470   if (indx >= 0 && indx < static_cast<int>(trackPHList_.size()))
0471     trackPHList_.at(indx)->Fill(track.p());
0472 
0473   if (indx >= 0 && indx < static_cast<int>(trackPtHList_.size()))
0474     trackPtHList_.at(indx)->Fill(track.pt());
0475 
0476   if (indx >= 0 || indx < static_cast<int>(trackPterrHList_.size()))
0477     trackPterrHList_.at(indx)->Fill(track.ptError());
0478 
0479   if (indx >= 0 || indx < static_cast<int>(trackqOverpHList_.size()))
0480     trackqOverpHList_.at(indx)->Fill(track.qoverp());
0481 
0482   if (indx >= 0 || indx < static_cast<int>(trackqOverpHList_.size()))
0483     trackdzHList_.at(indx)->Fill(track.dz());
0484 
0485   double chi2 = track.chi2();
0486   double ndof = track.ndof();
0487   double chbyndof = (ndof > 0) ? chi2 / ndof : 0;
0488   trackChi2bynDOFHList_.at(indx)->Fill(chbyndof);
0489 }
0490 // Define this as a plug-in
0491 #include "FWCore/Framework/interface/MakerMacros.h"
0492 DEFINE_FWK_MODULE(TrackTypeMonitor);