Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:18:42

0001 #include "HLTPixelIsolTrackFilter.h"
0002 
0003 #include "DataFormats/Common/interface/Handle.h"
0004 #include "DataFormats/Common/interface/RefToBase.h"
0005 
0006 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
0007 
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0011 
0012 HLTPixelIsolTrackFilter::HLTPixelIsolTrackFilter(const edm::ParameterSet& iConfig) : HLTFilter(iConfig) {
0013   candTag_ = iConfig.getParameter<edm::InputTag>("candTag");
0014   hltGTseedlabel_ = iConfig.getParameter<edm::InputTag>("L1GTSeedLabel");
0015   minDeltaPtL1Jet_ = iConfig.getParameter<double>("MinDeltaPtL1Jet");
0016   minpttrack_ = iConfig.getParameter<double>("MinPtTrack");
0017   maxptnearby_ = iConfig.getParameter<double>("MaxPtNearby");
0018   maxetatrack_ = iConfig.getParameter<double>("MaxEtaTrack");
0019   minetatrack_ = iConfig.getParameter<double>("MinEtaTrack");
0020   filterE_ = iConfig.getParameter<bool>("filterTrackEnergy");
0021   minEnergy_ = iConfig.getParameter<double>("MinEnergyTrack");
0022   nMaxTrackCandidates_ = iConfig.getParameter<int>("NMaxTrackCandidates");
0023   dropMultiL2Event_ = iConfig.getParameter<bool>("DropMultiL2Event");
0024   candToken_ = consumes<reco::IsolatedPixelTrackCandidateCollection>(candTag_);
0025   hltGTseedToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(hltGTseedlabel_);
0026 }
0027 
0028 HLTPixelIsolTrackFilter::~HLTPixelIsolTrackFilter() = default;
0029 
0030 void HLTPixelIsolTrackFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0031   edm::ParameterSetDescription desc;
0032   makeHLTFilterDescription(desc);
0033   desc.add<edm::InputTag>("candTag", edm::InputTag("isolPixelTrackProd"));
0034   desc.add<edm::InputTag>("L1GTSeedLabel", edm::InputTag("hltL1sIsoTrack"));
0035   desc.add<double>("MaxPtNearby", 2.0);
0036   desc.add<double>("MinEnergyTrack", 12.0);
0037   desc.add<double>("MinPtTrack", 3.5);
0038   desc.add<double>("MaxEtaTrack", 1.15);
0039   desc.add<double>("MinEtaTrack", 0.0);
0040   desc.add<double>("MinDeltaPtL1Jet", -40000.0);
0041   desc.add<bool>("filterTrackEnergy", true);
0042   desc.add<int>("NMaxTrackCandidates", 10);
0043   desc.add<bool>("DropMultiL2Event", false);
0044   descriptions.add("hltPixelIsolTrackFilter", desc);
0045 }
0046 
0047 bool HLTPixelIsolTrackFilter::hltFilter(edm::Event& iEvent,
0048                                         const edm::EventSetup& iSetup,
0049                                         trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0050   if (saveTags())
0051     filterproduct.addCollectionTag(candTag_);
0052 
0053   // Ref to Candidate object to be recorded in filter object
0054   edm::Ref<reco::IsolatedPixelTrackCandidateCollection> candref;
0055 
0056   // get hold of filtered candidates
0057   edm::Handle<reco::IsolatedPixelTrackCandidateCollection> recotrackcands;
0058   iEvent.getByToken(candToken_, recotrackcands);
0059 
0060   //Filtering
0061 
0062   //find leading L1 jet:
0063   double ptTriggered = -10;
0064 
0065   edm::Handle<trigger::TriggerFilterObjectWithRefs> l1trigobj;
0066   iEvent.getByToken(hltGTseedToken_, l1trigobj);
0067 
0068   std::vector<edm::Ref<l1extra::L1JetParticleCollection> > l1tauobjref;
0069   std::vector<edm::Ref<l1extra::L1JetParticleCollection> > l1jetobjref;
0070   std::vector<edm::Ref<l1extra::L1JetParticleCollection> > l1forjetobjref;
0071 
0072   l1trigobj->getObjects(trigger::TriggerL1TauJet, l1tauobjref);
0073   l1trigobj->getObjects(trigger::TriggerL1CenJet, l1jetobjref);
0074   l1trigobj->getObjects(trigger::TriggerL1ForJet, l1forjetobjref);
0075 
0076   for (auto& p : l1tauobjref)
0077     if (p->pt() > ptTriggered)
0078       ptTriggered = p->pt();
0079   for (auto& p : l1jetobjref)
0080     if (p->pt() > ptTriggered)
0081       ptTriggered = p->pt();
0082   for (auto& p : l1forjetobjref)
0083     if (p->pt() > ptTriggered)
0084       ptTriggered = p->pt();
0085 
0086   int n = 0;
0087   for (unsigned int i = 0; i < recotrackcands->size(); i++) {
0088     candref = edm::Ref<reco::IsolatedPixelTrackCandidateCollection>(recotrackcands, i);
0089 
0090     // cut on deltaPT
0091     if (ptTriggered - candref->pt() < minDeltaPtL1Jet_)
0092       continue;
0093 
0094     // select on transverse momentum
0095     if (!filterE_ && (candref->maxPtPxl() < maxptnearby_) && (candref->pt() > minpttrack_) &&
0096         fabs(candref->track()->eta()) < maxetatrack_ && fabs(candref->track()->eta()) > minetatrack_) {
0097       filterproduct.addObject(trigger::TriggerTrack, candref);
0098       n++;
0099       LogDebug("IsoTrk") << "PixelIsolP:Candidate[" << n << "] pt|eta|phi " << candref->pt() << "|" << candref->eta()
0100                          << "|" << candref->phi() << "\n";
0101     }
0102 
0103     // select on momentum
0104     if (filterE_) {
0105       if ((candref->maxPtPxl() < maxptnearby_) && ((candref->pt()) * cosh(candref->track()->eta()) > minEnergy_) &&
0106           fabs(candref->track()->eta()) < maxetatrack_ && fabs(candref->track()->eta()) > minetatrack_) {
0107         filterproduct.addObject(trigger::TriggerTrack, candref);
0108         n++;
0109         LogDebug("IsoTrk") << "PixelIsolE:Candidate[" << n << "] pt|eta|phi " << candref->pt() << "|" << candref->eta()
0110                            << "|" << candref->phi() << "\n";
0111       }
0112     }
0113 
0114     // stop looping over tracks if max number is reached
0115     if (!dropMultiL2Event_ && n >= nMaxTrackCandidates_)
0116       break;
0117 
0118   }  // loop over tracks
0119 
0120   bool accept(n > 0);
0121 
0122   if (dropMultiL2Event_ && n > nMaxTrackCandidates_)
0123     accept = false;
0124 
0125   return accept;
0126 }
0127 
0128 // declare this class as a framework plugin
0129 #include "FWCore/Framework/interface/MakerMacros.h"
0130 DEFINE_FWK_MODULE(HLTPixelIsolTrackFilter);