Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:27:44

0001 #include "RecoTauTag/HLTProducers/interface/TauJetSelectorForHLTTrackSeeding.h"
0002 
0003 #include "DataFormats/Math/interface/deltaPhi.h"
0004 
0005 TauJetSelectorForHLTTrackSeeding::TauJetSelectorForHLTTrackSeeding(const edm::ParameterSet& iConfig)
0006     : ptMinCaloJet_(iConfig.getParameter<double>("ptMinCaloJet")),
0007       etaMinCaloJet_(iConfig.getParameter<double>("etaMinCaloJet")),
0008       etaMaxCaloJet_(iConfig.getParameter<double>("etaMaxCaloJet")),
0009       tauConeSize_(iConfig.getParameter<double>("tauConeSize")),
0010       isolationConeSize_(iConfig.getParameter<double>("isolationConeSize")),
0011       fractionMinCaloInTauCone_(iConfig.getParameter<double>("fractionMinCaloInTauCone")),
0012       fractionMaxChargedPUInCaloCone_(iConfig.getParameter<double>("fractionMaxChargedPUInCaloCone")),
0013       ptTrkMaxInCaloCone_(iConfig.getParameter<double>("ptTrkMaxInCaloCone")),
0014       nTrkMaxInCaloCone_(iConfig.getParameter<int>("nTrkMaxInCaloCone")) {
0015   //now do what ever initialization is needed
0016   inputTrackJetToken_ = consumes<reco::TrackJetCollection>(iConfig.getParameter<edm::InputTag>("inputTrackJetTag"));
0017   inputCaloJetToken_ = consumes<reco::CaloJetCollection>(iConfig.getParameter<edm::InputTag>("inputCaloJetTag"));
0018   inputTrackToken_ = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("inputTrackTag"));
0019 
0020   produces<reco::TrackJetCollection>();
0021 }
0022 
0023 TauJetSelectorForHLTTrackSeeding::~TauJetSelectorForHLTTrackSeeding() {
0024   // do anything here that needs to be done at desctruction time
0025   // (e.g. close files, deallocate resources etc.)
0026 }
0027 
0028 //
0029 // member functions
0030 //
0031 
0032 // ------------ method called on each new Event  ------------
0033 void TauJetSelectorForHLTTrackSeeding::produce(edm::StreamID iStreamID,
0034                                                edm::Event& iEvent,
0035                                                const edm::EventSetup& iSetup) const {
0036   std::unique_ptr<reco::TrackJetCollection> augmentedTrackJets(new reco::TrackJetCollection);
0037 
0038   edm::Handle<reco::TrackJetCollection> trackjets;
0039   iEvent.getByToken(inputTrackJetToken_, trackjets);
0040 
0041   for (reco::TrackJetCollection::const_iterator trackjet = trackjets->begin(); trackjet != trackjets->end();
0042        trackjet++) {
0043     augmentedTrackJets->push_back(*trackjet);
0044   }
0045 
0046   edm::Handle<reco::TrackCollection> tracks;
0047   iEvent.getByToken(inputTrackToken_, tracks);
0048 
0049   edm::Handle<reco::CaloJetCollection> calojets;
0050   iEvent.getByToken(inputCaloJetToken_, calojets);
0051 
0052   const double tauConeSize2 = tauConeSize_ * tauConeSize_;
0053   const double isolationConeSize2 = isolationConeSize_ * isolationConeSize_;
0054 
0055   for (reco::CaloJetCollection::const_iterator calojet = calojets->begin(); calojet != calojets->end(); calojet++) {
0056     if (calojet->pt() < ptMinCaloJet_)
0057       continue;
0058     double etaJet = calojet->eta();
0059     double phiJet = calojet->phi();
0060     if (etaJet < etaMinCaloJet_)
0061       continue;
0062     if (etaJet > etaMaxCaloJet_)
0063       continue;
0064 
0065     std::vector<CaloTowerPtr> const& theTowers = calojet->getCaloConstituents();
0066     double ptIn = 0.;
0067     double ptOut = 0.;
0068     for (unsigned int itwr = 0; itwr < theTowers.size(); ++itwr) {
0069       double etaTwr = theTowers[itwr]->eta() - etaJet;
0070       double phiTwr = deltaPhi(theTowers[itwr]->phi(), phiJet);
0071       double deltaR2 = etaTwr * etaTwr + phiTwr * phiTwr;
0072       //std::cout << "Tower eta/phi/et : " << etaTwr << " " << phiTwr << " " << theTowers[itwr]->pt() << std::endl;
0073       if (deltaR2 < tauConeSize2) {
0074         ptIn += theTowers[itwr]->pt();
0075       } else if (deltaR2 < isolationConeSize2) {
0076         ptOut += theTowers[itwr]->pt();
0077       }
0078     }
0079     double ptTot = ptIn + ptOut;
0080     double fracIn = ptIn / ptTot;
0081 
0082     // We are looking for isolated tracks
0083     if (fracIn < fractionMinCaloInTauCone_)
0084       continue;
0085 
0086     int ntrk = 0;
0087     double ptTrk = 0.;
0088 
0089     for (reco::TrackJetCollection::const_iterator trackjet = trackjets->begin(); trackjet != trackjets->end();
0090          trackjet++) {
0091       for (unsigned itr = 0; itr < trackjet->numberOfTracks(); ++itr) {
0092         edm::Ptr<reco::Track> track = trackjet->track(itr);
0093         double trackEta = track->eta() - etaJet;
0094         double trackPhi = deltaPhi(track->phi(), phiJet);
0095         double deltaR2 = trackEta * trackEta + trackPhi * trackPhi;
0096         if (deltaR2 < isolationConeSize2) {
0097           ntrk++;
0098           ptTrk += track->pt();
0099         }
0100       }
0101     }
0102     // We are looking for calojets without signal tracks already in
0103     if (ntrk > nTrkMaxInCaloCone_)
0104       continue;
0105     if (ptTrk > ptTrkMaxInCaloCone_)
0106       continue;
0107 
0108     double ptTrk2 = 0.;
0109 
0110     for (reco::TrackCollection::const_iterator track = tracks->begin(); track != tracks->end(); track++) {
0111       double trackEta = track->eta() - etaJet;
0112       double trackPhi = deltaPhi(track->phi(), phiJet);
0113       double deltaR2 = trackEta * trackEta + trackPhi * trackPhi;
0114       if (deltaR2 < isolationConeSize2) {
0115         ptTrk2 += track->pt();
0116       }
0117     }
0118     // We are looking for signal jets, not PU jets
0119     double fractionChargedPU = ptTrk2 / calojet->pt();
0120     if (fractionChargedPU > fractionMaxChargedPUInCaloCone_)
0121       continue;
0122     /*
0123      std::cout << "Calo Jet " << calojet->pt() << " " << calojet->eta() 
0124            << " " << ptIn << " " << ptOut << " " << fracIn
0125            << " " << ptTrk << " " << ntrk 
0126            << " " << fractionChargedPU
0127            << std::endl;
0128      */
0129     math::XYZTLorentzVector p4(calojet->p4());
0130     math::XYZPoint vertex(calojet->vertex());
0131     augmentedTrackJets->push_back(reco::TrackJet(p4, vertex));
0132   }
0133 
0134   iEvent.put(std::move(augmentedTrackJets));
0135 }
0136 
0137 // ------------ method called once each job just before starting event loop  ------------
0138 void TauJetSelectorForHLTTrackSeeding::beginJob() {}
0139 
0140 // ------------ method called once each job just after ending the event loop  ------------
0141 void TauJetSelectorForHLTTrackSeeding::endJob() {}
0142 
0143 void TauJetSelectorForHLTTrackSeeding::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0144   edm::ParameterSetDescription desc;
0145   desc.add<edm::InputTag>("inputTrackJetTag", edm::InputTag("hltAntiKT5TrackJetsIter0"))
0146       ->setComment("Oryginal TrackJet collection");
0147   desc.add<edm::InputTag>("inputCaloJetTag", edm::InputTag("hltAntiKT5CaloJetsPFEt5"))
0148       ->setComment("CaloJet collection");
0149   desc.add<edm::InputTag>("inputTrackTag", edm::InputTag("hltPFlowTrackSelectionHighPurity"))
0150       ->setComment("Track collection for track isolation");
0151   desc.add<double>("ptMinCaloJet", 5.0)->setComment("Min. pT of CaloJet");
0152   desc.add<double>("etaMinCaloJet", -2.7)->setComment("Min. eta of CaloJet");
0153   desc.add<double>("etaMaxCaloJet", +2.7)->setComment("Max. eta of CaloJet");
0154   desc.add<double>("tauConeSize", 0.2)->setComment("Radius of tau signal cone");
0155   desc.add<double>("isolationConeSize", 0.5)->setComment("Radius of isolation cone");
0156   desc.add<double>("fractionMinCaloInTauCone", 0.8)->setComment("Min. fraction of calo energy in tau signal cone");
0157   desc.add<double>("fractionMaxChargedPUInCaloCone", 0.2)
0158       ->setComment("Max. fraction of PU charged energy Sum(Pt_trks)/Pt_jet in signal cone");
0159   desc.add<double>("ptTrkMaxInCaloCone", 1.0)->setComment("Max. sum of pT of tracks in isolation cone");
0160   desc.add<int>("nTrkMaxInCaloCone", 0)->setComment("Max. no. of tracks in isolation cone");
0161   descriptions.setComment(
0162       "This module produces a collection of TrackJets to seed iterative tracking.\nIt is done by enriching the input "
0163       "TrackJet collection with CaloJets passing isolation critera - tau candidates");
0164   descriptions.add("tauJetSelectorForHLTTrackSeeding", desc);
0165 }
0166 
0167 #include "FWCore/Framework/interface/MakerMacros.h"
0168 DEFINE_FWK_MODULE(TauJetSelectorForHLTTrackSeeding);