Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-10-08 05:12:15

0001 #include "Validation/RecoB/plugins/BDHadronTrackMonitoringAnalyzer.h"
0002 
0003 using namespace reco;
0004 using namespace edm;
0005 using namespace std;
0006 
0007 const reco::TrackBaseRef toTrackRef(const edm::Ptr<reco::Candidate> &cnd) {
0008   const reco::PFCandidate *pfcand = dynamic_cast<const reco::PFCandidate *>(cnd.get());
0009 
0010   if ((std::abs(pfcand->pdgId()) == 11 || pfcand->pdgId() == 22) && pfcand->gsfTrackRef().isNonnull() &&
0011       pfcand->gsfTrackRef().isAvailable())
0012     return reco::TrackBaseRef(pfcand->gsfTrackRef());
0013   else if (pfcand->trackRef().isNonnull() && pfcand->trackRef().isAvailable())
0014     return reco::TrackBaseRef(pfcand->trackRef());
0015   else
0016     return reco::TrackBaseRef();
0017 }
0018 
0019 // ---------- Static member declaration -----------
0020 const std::vector<std::string> BDHadronTrackMonitoringAnalyzer::TrkHistCat = {
0021     "BCWeakDecay", "BWeakDecay", "CWeakDecay", "PU", "Other", "Fake"};
0022 
0023 // ---------- Constructor -----------
0024 BDHadronTrackMonitoringAnalyzer::BDHadronTrackMonitoringAnalyzer(const edm::ParameterSet &pSet)
0025     : distJetAxis_(pSet.getParameter<double>("distJetAxisCut")),
0026       decayLength_(pSet.getParameter<double>("decayLengthCut")),
0027       minJetPt_(pSet.getParameter<double>("minJetPt")),
0028       maxJetEta_(pSet.getParameter<double>("maxJetEta")),
0029       ipTagInfos_(pSet.getParameter<std::string>("ipTagInfos")),
0030       PatJetSrc_(pSet.getParameter<InputTag>("PatJetSource")),
0031       TrackSrc_(pSet.getParameter<InputTag>("TrackSource")),
0032       PVSrc_(pSet.getParameter<InputTag>("PrimaryVertexSource")),
0033       ClusterTPMapSrc_(pSet.getParameter<InputTag>("clusterTPMap")),
0034       classifier_(pSet, consumesCollector())
0035 
0036 {
0037   PatJetCollectionTag_ = consumes<pat::JetCollection>(PatJetSrc_);
0038   TrackCollectionTag_ = consumes<reco::TrackCollection>(TrackSrc_);
0039   PrimaryVertexColl_ = consumes<reco::VertexCollection>(PVSrc_);
0040   clusterTPMapToken_ = consumes<ClusterTPAssociation>(ClusterTPMapSrc_);
0041   ttrackToken_ = esConsumes(edm::ESInputTag("", "TransientTrackBuilder"));
0042   // TrkHistCat = {"BCWeakDecay", "BWeakDecay", "CWeakDecay", "PU", "Other",
0043   // "Fake"};
0044 }
0045 
0046 // ---------- BookHistograms -----------
0047 
0048 void BDHadronTrackMonitoringAnalyzer::bookHistograms(DQMStore::IBooker &ibook,
0049                                                      edm::Run const &run,
0050                                                      edm::EventSetup const &es) {
0051   ibook.setCurrentFolder("BDHadronTracks/JetContent");
0052   //
0053   // Book all histograms.
0054   //
0055   RecoBTag::setTDRStyle();
0056 
0057   nTrkAll_bjet = ibook.book1D(
0058       "nTrkAll_bjet", "Number of selected tracks in b jets;number of selected tracks;jets", 16, -0.5, 15.5);
0059 
0060   nTrkAll_cjet = ibook.book1D(
0061       "nTrkAll_cjet", "Number of selected tracks in c jets;number of selected tracks;jets", 16, -0.5, 15.5);
0062 
0063   nTrkAll_dusgjet = ibook.book1D(
0064       "nTrkAll_dusgjet", "Number of selected tracks in dusg jets;number of selected tracks;jets", 16, -0.5, 15.5);
0065 
0066   // Loop over different Track History Categories
0067   for (unsigned int i = 0; i < TrkHistCat.size(); i++) {
0068     ibook.setCurrentFolder("BDHadronTracks/JetContent");
0069     // b jets
0070     nTrk_bjet[i] = ibook.book1D("nTrk_bjet_" + TrkHistCat[i],
0071                                 "Number of selected tracks in b jets (" + TrkHistCat[i] +
0072                                     ");number of selected tracks (" + TrkHistCat[i] + ");jets",
0073                                 16,
0074                                 -0.5,
0075                                 15.5);
0076 
0077     // c jets
0078     nTrk_cjet[i] = ibook.book1D("nTrk_cjet_" + TrkHistCat[i],
0079                                 "Number of selected tracks in c jets (" + TrkHistCat[i] +
0080                                     ");number of selected tracks (" + TrkHistCat[i] + ");jets",
0081                                 16,
0082                                 -0.5,
0083                                 15.5);
0084 
0085     // dusg jets
0086     nTrk_dusgjet[i] = ibook.book1D("nTrk_dusgjet_" + TrkHistCat[i],
0087                                    "Number of selected tracks in dusg jets (" + TrkHistCat[i] +
0088                                        ");number of selected tracks (" + TrkHistCat[i] + ");jets",
0089                                    16,
0090                                    -0.5,
0091                                    15.5);
0092 
0093     ibook.setCurrentFolder("BDHadronTracks/TrackInfo");
0094     // track properties for all flavours combined
0095     TrkPt_alljets[i] = ibook.book1D("TrkPt_" + TrkHistCat[i],
0096                                     "Track pT (" + TrkHistCat[i] + ");track p_{T} (" + TrkHistCat[i] + ");tracks",
0097                                     30,
0098                                     0,
0099                                     100);
0100     TrkEta_alljets[i] = ibook.book1D("TrkEta_" + TrkHistCat[i],
0101                                      "Track #eta (" + TrkHistCat[i] + ");track #eta (" + TrkHistCat[i] + ");tracks",
0102                                      30,
0103                                      -2.5,
0104                                      2.5);
0105     TrkPhi_alljets[i] = ibook.book1D("TrkPhi_" + TrkHistCat[i],
0106                                      "Track #phi (" + TrkHistCat[i] + ");track #phi (" + TrkHistCat[i] + ");tracks",
0107                                      30,
0108                                      -3.15,
0109                                      3.15);
0110     TrkDxy_alljets[i] = ibook.book1D("TrkDxy_" + TrkHistCat[i],
0111                                      "Track dxy (" + TrkHistCat[i] + ");track dxy (" + TrkHistCat[i] + ");tracks",
0112                                      30,
0113                                      -0.1,
0114                                      0.1);
0115     TrkDz_alljets[i] = ibook.book1D("TrkDz_" + TrkHistCat[i],
0116                                     "Track dz (" + TrkHistCat[i] + ");track dz (" + TrkHistCat[i] + ");tracks",
0117                                     30,
0118                                     -0.1,
0119                                     0.1);
0120     TrkHitAll_alljets[i] = ibook.book1D(
0121         "TrkHitAll_" + TrkHistCat[i],
0122         "Number of tracker hits (" + TrkHistCat[i] + ");track number of all hits (" + TrkHistCat[i] + ");tracks",
0123         31,
0124         -0.5,
0125         30.5);
0126     TrkHitStrip_alljets[i] = ibook.book1D(
0127         "TrkHitStrip_" + TrkHistCat[i],
0128         "Number of strip hits (" + TrkHistCat[i] + ");track number of strip hits (" + TrkHistCat[i] + ");tracks",
0129         31,
0130         -0.5,
0131         30.5);
0132     TrkHitPixel_alljets[i] = ibook.book1D(
0133         "TrkHitPixel_" + TrkHistCat[i],
0134         "Number of pixel hits (" + TrkHistCat[i] + ");track number of pixel hits (" + TrkHistCat[i] + ");tracks",
0135         9,
0136         -0.5,
0137         8.5);
0138 
0139     ibook.setCurrentFolder("BDHadronTracks/TrackTruthInfo");
0140     if (i < 5) {  // Fakes (i == 5) have no truth by definition!
0141       TrkTruthPt_alljets[i] =
0142           ibook.book1D("TrkTruthPt_" + TrkHistCat[i],
0143                        "Track pT (" + TrkHistCat[i] + " Truth);track p_{T} (" + TrkHistCat[i] + " Truth);tracks",
0144                        30,
0145                        0,
0146                        100);
0147       TrkTruthEta_alljets[i] =
0148           ibook.book1D("TrkTruthEta_" + TrkHistCat[i],
0149                        "Track #eta (" + TrkHistCat[i] + " Truth);track #eta (" + TrkHistCat[i] + " Truth);tracks",
0150                        30,
0151                        -2.5,
0152                        2.5);
0153       TrkTruthPhi_alljets[i] =
0154           ibook.book1D("TrkTruthPhi_" + TrkHistCat[i],
0155                        "Track #phi (" + TrkHistCat[i] + " Truth);track #phi (" + TrkHistCat[i] + " Truth);tracks",
0156                        30,
0157                        -3.15,
0158                        3.15);
0159       TrkTruthDxy_alljets[i] =
0160           ibook.book1D("TrkTruthDxy_" + TrkHistCat[i],
0161                        "Track dxy (" + TrkHistCat[i] + " Truth);track dxy (" + TrkHistCat[i] + " Truth);tracks",
0162                        30,
0163                        -0.1,
0164                        0.1);
0165       TrkTruthDz_alljets[i] =
0166           ibook.book1D("TrkTruthDz_" + TrkHistCat[i],
0167                        "Track dz (" + TrkHistCat[i] + " Truth);track dz (" + TrkHistCat[i] + " Truth);tracks",
0168                        30,
0169                        -0.1,
0170                        0.1);
0171       TrkTruthHitAll_alljets[i] =
0172           ibook.book1D("TrkTruthHitAll_" + TrkHistCat[i],
0173                        "Number of tracker hits (" + TrkHistCat[i] + " Truth);track number of all hits (" +
0174                            TrkHistCat[i] + " Truth);tracks",
0175                        31,
0176                        -0.5,
0177                        30.5);
0178       TrkTruthHitStrip_alljets[i] =
0179           ibook.book1D("TrkTruthHitStrip_" + TrkHistCat[i],
0180                        "Number of strip hits (" + TrkHistCat[i] + " Truth);track number of strip hits (" +
0181                            TrkHistCat[i] + " Truth);tracks",
0182                        31,
0183                        -0.5,
0184                        30.5);
0185       TrkTruthHitPixel_alljets[i] =
0186           ibook.book1D("TrkTruthHitPixel_" + TrkHistCat[i],
0187                        "Number of pixel hits (" + TrkHistCat[i] + " Truth);track number of pixel hits (" +
0188                            TrkHistCat[i] + " Truth);tracks",
0189                        9,
0190                        -0.5,
0191                        8.5);
0192     }
0193   }
0194 }
0195 
0196 // ---------- Destructor -----------
0197 
0198 BDHadronTrackMonitoringAnalyzer::~BDHadronTrackMonitoringAnalyzer() {}
0199 
0200 // ---------- Analyze -----------
0201 // This is needed to get a TrackingParticle --> Cluster match (instead of
0202 // Cluster-->TP)
0203 using P = std::pair<OmniClusterRef, TrackingParticleRef>;
0204 bool compare(const P &i, const P &j) { return i.second.index() > j.second.index(); }
0205 
0206 void BDHadronTrackMonitoringAnalyzer::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0207   edm::Handle<pat::JetCollection> patJetsColl;
0208   iEvent.getByToken(PatJetCollectionTag_, patJetsColl);
0209 
0210   edm::Handle<reco::TrackCollection> tracksHandle;
0211   iEvent.getByToken(TrackCollectionTag_, tracksHandle);
0212 
0213   edm::Handle<ClusterTPAssociation> pCluster2TPListH;
0214   iEvent.getByToken(clusterTPMapToken_, pCluster2TPListH);
0215   const ClusterTPAssociation &clusterToTPMap = *pCluster2TPListH;
0216 
0217   //edm::ESHandle<TransientTrackBuilder>
0218   const auto &trackBuilder = iSetup.getHandle(ttrackToken_);
0219 
0220   classifier_.newEvent(iEvent, iSetup);
0221 
0222   // -----Primary Vertex-----
0223   const reco::Vertex *pv;
0224 
0225   edm::Handle<reco::VertexCollection> primaryVertex;
0226   iEvent.getByToken(PrimaryVertexColl_, primaryVertex);
0227 
0228   bool pvFound = (!primaryVertex->empty());
0229   if (pvFound) {
0230     pv = &(*primaryVertex->begin());
0231   } else {
0232     reco::Vertex::Error e;
0233     e(0, 0) = 0.0015 * 0.0015;
0234     e(1, 1) = 0.0015 * 0.0015;
0235     e(2, 2) = 15. * 15.;
0236     reco::Vertex::Point p(0, 0, 0);
0237     pv = new reco::Vertex(p, e, 1, 1, 1);
0238   }
0239   // -----------------------
0240 
0241   // -------- Loop Over Jets ----------
0242   for (pat::JetCollection::const_iterator jet = patJetsColl->begin(); jet != patJetsColl->end(); ++jet) {
0243     if (jet->pt() < minJetPt_ || std::fabs(jet->eta()) > maxJetEta_)
0244       continue;
0245 
0246     unsigned int flav = abs(jet->hadronFlavour());
0247 
0248     // std::cout << "patJet collection has pfImpactParameterTagInfo?: " <<
0249     // jet->hasTagInfo("pfImpactParameter") << std::endl;
0250     const CandIPTagInfo *trackIpTagInfo = jet->tagInfoCandIP(ipTagInfos_);
0251     const std::vector<edm::Ptr<reco::Candidate>> &selectedTracks(trackIpTagInfo->selectedTracks());
0252 
0253     unsigned int nseltracks = 0;
0254     std::vector<int> nseltracksCat(TrkHistCat.size(),
0255                                    0);  // following the order of TrkHistCat
0256 
0257     unsigned int nTrackSize = selectedTracks.size();  // number of tracks from IPInfos to loop over
0258     // -------- Loop Over (selected) Tracks ----------
0259     for (unsigned int itt = 0; itt < nTrackSize; ++itt) {
0260       const TrackBaseRef ptrackRef = toTrackRef(selectedTracks[itt]);
0261       const reco::Track *ptrackPtr = reco::btag::toTrack(ptrackRef);
0262       const reco::Track &ptrack = *ptrackPtr;
0263 
0264       reco::TransientTrack transientTrack = trackBuilder->build(ptrackPtr);
0265       GlobalVector direction(jet->px(), jet->py(), jet->pz());
0266 
0267       Double_t distJetAxis = IPTools::jetTrackDistance(transientTrack, direction, *pv).second.value();
0268 
0269       Double_t decayLength = 999;
0270       TrajectoryStateOnSurface closest =
0271           IPTools::closestApproachToJet(transientTrack.impactPointState(), *pv, direction, transientTrack.field());
0272       if (closest.isValid())
0273         decayLength = (closest.globalPosition() - RecoVertex::convertPos(pv->position())).mag();
0274       else
0275         decayLength = 999;
0276 
0277       // extra cut ons the tracks
0278       if (std::fabs(distJetAxis) > distJetAxis_ || decayLength > decayLength_) {
0279         continue;
0280       }
0281       nseltracks += 1;  // if it passed these cuts, nselectedtracks +1
0282 
0283       TrackCategories::Flags theFlag = classifier_.evaluate(toTrackRef(selectedTracks[itt])).flags();
0284 
0285       double TrkPt = ptrack.pt();
0286       double TrkEta = ptrack.eta();
0287       double TrkPhi = ptrack.phi();
0288       double TrkDxy = ptrack.dxy(pv->position());
0289       double TrkDz = ptrack.dz(pv->position());
0290       int TrknHitAll = ptrack.numberOfValidHits();
0291       int TrknHitPixel = ptrack.hitPattern().numberOfValidPixelHits();
0292       int TrknHitStrip = ptrack.hitPattern().numberOfValidStripHits();
0293 
0294       double TrkTruthPt = -99;
0295       double TrkTruthEta = -99;
0296       double TrkTruthPhi = -99;
0297       double TrkTruthDxy = -1;
0298       double TrkTruthDz = -1;
0299       int TrkTruthnHitAll = -1;
0300       int TrkTruthnHitPixel = -1;
0301       int TrkTruthnHitStrip = -1;
0302 
0303       // Get corresponding Trackingparticle
0304       std::pair<TrackingParticleRef, double> res = classifier_.history().getMatchedTrackingParticle();
0305       TrackingParticleRef tpr = res.first;
0306       double quality_tpr = res.second;
0307 
0308       // Match TP to hit-cluster (re-ordering according to TP rather than
0309       // clusters and look for equal_range of a given tpr)
0310       auto clusterTPmap = clusterToTPMap.map();
0311       std::sort(clusterTPmap.begin(), clusterTPmap.end(), compare);
0312       auto clusterRange =
0313           std::equal_range(clusterTPmap.begin(), clusterTPmap.end(), std::make_pair(OmniClusterRef(), tpr), compare);
0314       if (quality_tpr != 0) {
0315         TrkTruthPt = tpr->pt();
0316         TrkTruthEta = tpr->eta();
0317         TrkTruthPhi = tpr->phi();
0318 
0319         const TrackingParticle::Point &vertex_pv = pv->position();
0320         TrackingParticle::Point vertex_tpr = tpr->vertex();
0321         TrackingParticle::Vector momentum_tpr = tpr->momentum();
0322         TrkTruthDxy = (-(vertex_tpr.x() - vertex_pv.x()) * momentum_tpr.y() +
0323                        (vertex_tpr.y() - vertex_pv.y()) * momentum_tpr.x()) /
0324                       tpr->pt();
0325         TrkTruthDz = (vertex_tpr.z() - vertex_pv.z()) - ((vertex_tpr.x() - vertex_pv.x()) * momentum_tpr.x() +
0326                                                          (vertex_tpr.y() - vertex_pv.y()) * momentum_tpr.y()) /
0327                                                             sqrt(momentum_tpr.perp2()) * momentum_tpr.z() /
0328                                                             sqrt(momentum_tpr.perp2());
0329 
0330         [[clang::suppress]] TrkTruthnHitAll = 0;
0331         TrkTruthnHitPixel = 0;
0332         TrkTruthnHitStrip = 0;
0333         if (clusterRange.first != clusterRange.second) {
0334           for (auto ip = clusterRange.first; ip != clusterRange.second; ++ip) {
0335             const OmniClusterRef &cluster = ip->first;
0336             if (cluster.isPixel() && cluster.isValid()) {
0337               TrkTruthnHitPixel += 1;
0338             }
0339             if (cluster.isStrip() && cluster.isValid()) {
0340               TrkTruthnHitStrip += 1;
0341             }
0342           }
0343         }
0344         TrkTruthnHitAll = TrkTruthnHitPixel + TrkTruthnHitStrip;
0345       }
0346 
0347       // ----------- Filling the correct histograms based on jet flavour and
0348       // Track history Category --------
0349 
0350       // BCWeakDecay
0351       if (theFlag[TrackCategories::SignalEvent] && theFlag[TrackCategories::BWeakDecay] &&
0352           theFlag[TrackCategories::CWeakDecay]) {
0353         nseltracksCat[BDHadronTrackMonitoringAnalyzer::BCWeakDecay] += 1;
0354         TrkPt_alljets[BDHadronTrackMonitoringAnalyzer::BCWeakDecay]->Fill(TrkPt);
0355         TrkEta_alljets[BDHadronTrackMonitoringAnalyzer::BCWeakDecay]->Fill(TrkEta);
0356         TrkPhi_alljets[BDHadronTrackMonitoringAnalyzer::BCWeakDecay]->Fill(TrkPhi);
0357         TrkDxy_alljets[BDHadronTrackMonitoringAnalyzer::BCWeakDecay]->Fill(TrkDxy);
0358         TrkDz_alljets[BDHadronTrackMonitoringAnalyzer::BCWeakDecay]->Fill(TrkDz);
0359         TrkHitAll_alljets[BDHadronTrackMonitoringAnalyzer::BCWeakDecay]->Fill(TrknHitAll);
0360         TrkHitPixel_alljets[BDHadronTrackMonitoringAnalyzer::BCWeakDecay]->Fill(TrknHitPixel);
0361         TrkHitStrip_alljets[BDHadronTrackMonitoringAnalyzer::BCWeakDecay]->Fill(TrknHitStrip);
0362         if (quality_tpr != 0) {
0363           TrkTruthPt_alljets[BDHadronTrackMonitoringAnalyzer::BCWeakDecay]->Fill(TrkTruthPt);
0364           TrkTruthEta_alljets[BDHadronTrackMonitoringAnalyzer::BCWeakDecay]->Fill(TrkTruthEta);
0365           TrkTruthPhi_alljets[BDHadronTrackMonitoringAnalyzer::BCWeakDecay]->Fill(TrkTruthPhi);
0366           TrkTruthDxy_alljets[BDHadronTrackMonitoringAnalyzer::BCWeakDecay]->Fill(TrkTruthDxy);
0367           TrkTruthDz_alljets[BDHadronTrackMonitoringAnalyzer::BCWeakDecay]->Fill(TrkTruthDz);
0368           TrkTruthHitAll_alljets[BDHadronTrackMonitoringAnalyzer::BCWeakDecay]->Fill(TrkTruthnHitAll);
0369           TrkTruthHitPixel_alljets[BDHadronTrackMonitoringAnalyzer::BCWeakDecay]->Fill(TrkTruthnHitPixel);
0370           TrkTruthHitStrip_alljets[BDHadronTrackMonitoringAnalyzer::BCWeakDecay]->Fill(TrkTruthnHitStrip);
0371         }
0372       }
0373       // BWeakDecay
0374       else if (theFlag[TrackCategories::SignalEvent] && theFlag[TrackCategories::BWeakDecay] &&
0375                !theFlag[TrackCategories::CWeakDecay]) {
0376         nseltracksCat[BDHadronTrackMonitoringAnalyzer::BWeakDecay] += 1;
0377         TrkPt_alljets[BDHadronTrackMonitoringAnalyzer::BWeakDecay]->Fill(TrkPt);
0378         TrkEta_alljets[BDHadronTrackMonitoringAnalyzer::BWeakDecay]->Fill(TrkEta);
0379         TrkPhi_alljets[BDHadronTrackMonitoringAnalyzer::BWeakDecay]->Fill(TrkPhi);
0380         TrkDxy_alljets[BDHadronTrackMonitoringAnalyzer::BWeakDecay]->Fill(TrkDxy);
0381         TrkDz_alljets[BDHadronTrackMonitoringAnalyzer::BWeakDecay]->Fill(TrkDz);
0382         TrkHitAll_alljets[BDHadronTrackMonitoringAnalyzer::BWeakDecay]->Fill(TrknHitAll);
0383         TrkHitPixel_alljets[BDHadronTrackMonitoringAnalyzer::BWeakDecay]->Fill(TrknHitPixel);
0384         TrkHitStrip_alljets[BDHadronTrackMonitoringAnalyzer::BWeakDecay]->Fill(TrknHitStrip);
0385         if (quality_tpr != 0) {
0386           TrkTruthPt_alljets[BDHadronTrackMonitoringAnalyzer::BWeakDecay]->Fill(TrkTruthPt);
0387           TrkTruthEta_alljets[BDHadronTrackMonitoringAnalyzer::BWeakDecay]->Fill(TrkTruthEta);
0388           TrkTruthPhi_alljets[BDHadronTrackMonitoringAnalyzer::BWeakDecay]->Fill(TrkTruthPhi);
0389           TrkTruthDxy_alljets[BDHadronTrackMonitoringAnalyzer::BWeakDecay]->Fill(TrkTruthDxy);
0390           TrkTruthDz_alljets[BDHadronTrackMonitoringAnalyzer::BWeakDecay]->Fill(TrkTruthDz);
0391           TrkTruthHitAll_alljets[BDHadronTrackMonitoringAnalyzer::BWeakDecay]->Fill(TrkTruthnHitAll);
0392           TrkTruthHitPixel_alljets[BDHadronTrackMonitoringAnalyzer::BWeakDecay]->Fill(TrkTruthnHitPixel);
0393           TrkTruthHitStrip_alljets[BDHadronTrackMonitoringAnalyzer::BWeakDecay]->Fill(TrkTruthnHitStrip);
0394         }
0395       }
0396       // CWeakDecay
0397       else if (theFlag[TrackCategories::SignalEvent] && !theFlag[TrackCategories::BWeakDecay] &&
0398                theFlag[TrackCategories::CWeakDecay]) {
0399         nseltracksCat[BDHadronTrackMonitoringAnalyzer::CWeakDecay] += 1;
0400         TrkPt_alljets[BDHadronTrackMonitoringAnalyzer::CWeakDecay]->Fill(TrkPt);
0401         TrkEta_alljets[BDHadronTrackMonitoringAnalyzer::CWeakDecay]->Fill(TrkEta);
0402         TrkPhi_alljets[BDHadronTrackMonitoringAnalyzer::CWeakDecay]->Fill(TrkPhi);
0403         TrkDxy_alljets[BDHadronTrackMonitoringAnalyzer::CWeakDecay]->Fill(TrkDxy);
0404         TrkDz_alljets[BDHadronTrackMonitoringAnalyzer::CWeakDecay]->Fill(TrkDz);
0405         TrkHitAll_alljets[BDHadronTrackMonitoringAnalyzer::CWeakDecay]->Fill(TrknHitAll);
0406         TrkHitPixel_alljets[BDHadronTrackMonitoringAnalyzer::CWeakDecay]->Fill(TrknHitPixel);
0407         TrkHitStrip_alljets[BDHadronTrackMonitoringAnalyzer::CWeakDecay]->Fill(TrknHitStrip);
0408         if (quality_tpr != 0) {
0409           TrkTruthPt_alljets[BDHadronTrackMonitoringAnalyzer::CWeakDecay]->Fill(TrkTruthPt);
0410           TrkTruthEta_alljets[BDHadronTrackMonitoringAnalyzer::CWeakDecay]->Fill(TrkTruthEta);
0411           TrkTruthPhi_alljets[BDHadronTrackMonitoringAnalyzer::CWeakDecay]->Fill(TrkTruthPhi);
0412           TrkTruthDxy_alljets[BDHadronTrackMonitoringAnalyzer::CWeakDecay]->Fill(TrkTruthDxy);
0413           TrkTruthDz_alljets[BDHadronTrackMonitoringAnalyzer::CWeakDecay]->Fill(TrkTruthDz);
0414           TrkTruthHitAll_alljets[BDHadronTrackMonitoringAnalyzer::CWeakDecay]->Fill(TrkTruthnHitAll);
0415           TrkTruthHitPixel_alljets[BDHadronTrackMonitoringAnalyzer::CWeakDecay]->Fill(TrkTruthnHitPixel);
0416           TrkTruthHitStrip_alljets[BDHadronTrackMonitoringAnalyzer::CWeakDecay]->Fill(TrkTruthnHitStrip);
0417         }
0418       }
0419       // PU
0420       else if (!theFlag[TrackCategories::SignalEvent] && !theFlag[TrackCategories::Fake]) {
0421         nseltracksCat[BDHadronTrackMonitoringAnalyzer::PU] += 1;
0422         TrkPt_alljets[BDHadronTrackMonitoringAnalyzer::PU]->Fill(TrkPt);
0423         TrkEta_alljets[BDHadronTrackMonitoringAnalyzer::PU]->Fill(TrkEta);
0424         TrkPhi_alljets[BDHadronTrackMonitoringAnalyzer::PU]->Fill(TrkPhi);
0425         TrkDxy_alljets[BDHadronTrackMonitoringAnalyzer::PU]->Fill(TrkDxy);
0426         TrkDz_alljets[BDHadronTrackMonitoringAnalyzer::PU]->Fill(TrkDz);
0427         TrkHitAll_alljets[BDHadronTrackMonitoringAnalyzer::PU]->Fill(TrknHitAll);
0428         TrkHitPixel_alljets[BDHadronTrackMonitoringAnalyzer::PU]->Fill(TrknHitPixel);
0429         TrkHitStrip_alljets[BDHadronTrackMonitoringAnalyzer::PU]->Fill(TrknHitStrip);
0430         if (quality_tpr != 0) {
0431           TrkTruthPt_alljets[BDHadronTrackMonitoringAnalyzer::PU]->Fill(TrkTruthPt);
0432           TrkTruthEta_alljets[BDHadronTrackMonitoringAnalyzer::PU]->Fill(TrkTruthEta);
0433           TrkTruthPhi_alljets[BDHadronTrackMonitoringAnalyzer::PU]->Fill(TrkTruthPhi);
0434           TrkTruthDxy_alljets[BDHadronTrackMonitoringAnalyzer::PU]->Fill(TrkTruthDxy);
0435           TrkTruthDz_alljets[BDHadronTrackMonitoringAnalyzer::PU]->Fill(TrkTruthDz);
0436           TrkTruthHitAll_alljets[BDHadronTrackMonitoringAnalyzer::PU]->Fill(TrkTruthnHitAll);
0437           TrkTruthHitPixel_alljets[BDHadronTrackMonitoringAnalyzer::PU]->Fill(TrkTruthnHitPixel);
0438           TrkTruthHitStrip_alljets[BDHadronTrackMonitoringAnalyzer::PU]->Fill(TrkTruthnHitStrip);
0439         }
0440       }
0441       // Other
0442       else if (theFlag[TrackCategories::SignalEvent] && !theFlag[TrackCategories::BWeakDecay] &&
0443                !theFlag[TrackCategories::CWeakDecay]) {
0444         nseltracksCat[BDHadronTrackMonitoringAnalyzer::Other] += 1;
0445         TrkPt_alljets[BDHadronTrackMonitoringAnalyzer::Other]->Fill(TrkPt);
0446         TrkEta_alljets[BDHadronTrackMonitoringAnalyzer::Other]->Fill(TrkEta);
0447         TrkPhi_alljets[BDHadronTrackMonitoringAnalyzer::Other]->Fill(TrkPhi);
0448         TrkDxy_alljets[BDHadronTrackMonitoringAnalyzer::Other]->Fill(TrkDxy);
0449         TrkDz_alljets[BDHadronTrackMonitoringAnalyzer::Other]->Fill(TrkDz);
0450         TrkHitAll_alljets[BDHadronTrackMonitoringAnalyzer::Other]->Fill(TrknHitAll);
0451         TrkHitPixel_alljets[BDHadronTrackMonitoringAnalyzer::Other]->Fill(TrknHitPixel);
0452         TrkHitStrip_alljets[BDHadronTrackMonitoringAnalyzer::Other]->Fill(TrknHitStrip);
0453         if (quality_tpr != 0) {
0454           TrkTruthPt_alljets[BDHadronTrackMonitoringAnalyzer::Other]->Fill(TrkTruthPt);
0455           TrkTruthEta_alljets[BDHadronTrackMonitoringAnalyzer::Other]->Fill(TrkTruthEta);
0456           TrkTruthPhi_alljets[BDHadronTrackMonitoringAnalyzer::Other]->Fill(TrkTruthPhi);
0457           TrkTruthDxy_alljets[BDHadronTrackMonitoringAnalyzer::Other]->Fill(TrkTruthDxy);
0458           TrkTruthDz_alljets[BDHadronTrackMonitoringAnalyzer::Other]->Fill(TrkTruthDz);
0459           TrkTruthHitAll_alljets[BDHadronTrackMonitoringAnalyzer::Other]->Fill(TrkTruthnHitAll);
0460           TrkTruthHitPixel_alljets[BDHadronTrackMonitoringAnalyzer::Other]->Fill(TrkTruthnHitPixel);
0461           TrkTruthHitStrip_alljets[BDHadronTrackMonitoringAnalyzer::Other]->Fill(TrkTruthnHitStrip);
0462         }
0463       }
0464       // Fake
0465       else if (!theFlag[TrackCategories::SignalEvent] && theFlag[TrackCategories::Fake]) {
0466         nseltracksCat[BDHadronTrackMonitoringAnalyzer::Fake] += 1;
0467         TrkPt_alljets[BDHadronTrackMonitoringAnalyzer::Fake]->Fill(TrkPt);
0468         TrkEta_alljets[BDHadronTrackMonitoringAnalyzer::Fake]->Fill(TrkEta);
0469         TrkPhi_alljets[BDHadronTrackMonitoringAnalyzer::Fake]->Fill(TrkPhi);
0470         TrkDxy_alljets[BDHadronTrackMonitoringAnalyzer::Fake]->Fill(TrkDxy);
0471         TrkDz_alljets[BDHadronTrackMonitoringAnalyzer::Fake]->Fill(TrkDz);
0472         TrkHitAll_alljets[BDHadronTrackMonitoringAnalyzer::Fake]->Fill(TrknHitAll);
0473         TrkHitPixel_alljets[BDHadronTrackMonitoringAnalyzer::Fake]->Fill(TrknHitPixel);
0474         TrkHitStrip_alljets[BDHadronTrackMonitoringAnalyzer::Fake]->Fill(TrknHitStrip);
0475         // NO TRUTH FOR FAKES!!!
0476       }
0477     }
0478     // -------- END Loop Over (selected) Tracks ----------
0479     // Still have to fill some jet-flavour specific variables
0480     if (flav == 5) {
0481       nTrkAll_bjet->Fill(nseltracks);
0482       for (unsigned int i = 0; i < TrkHistCat.size(); i++) {
0483         nTrk_bjet[i]->Fill(nseltracksCat[i]);
0484       }
0485     } else if (flav == 4) {
0486       nTrkAll_cjet->Fill(nseltracks);
0487       for (unsigned int i = 0; i < TrkHistCat.size(); i++) {
0488         nTrk_cjet[i]->Fill(nseltracksCat[i]);
0489       }
0490     } else {
0491       nTrkAll_dusgjet->Fill(nseltracks);
0492       for (unsigned int i = 0; i < TrkHistCat.size(); i++) {
0493         nTrk_dusgjet[i]->Fill(nseltracksCat[i]);
0494       }
0495     }
0496   }
0497   // -------- END Loop Over Jets ----------
0498 
0499   if (!pvFound) {
0500     delete pv;
0501   }
0502 }
0503 
0504 // define this as a plug-in
0505 DEFINE_FWK_MODULE(BDHadronTrackMonitoringAnalyzer);