CMS 3D CMS Logo

HLTPixelIsolTrackFilter.cc
Go to the documentation of this file.
2 
5 
7 
11 
13  candTag_ = iConfig.getParameter<edm::InputTag> ("candTag");
14  hltGTseedlabel_ = iConfig.getParameter<edm::InputTag> ("L1GTSeedLabel");
15  minDeltaPtL1Jet_ = iConfig.getParameter<double> ("MinDeltaPtL1Jet");
16  minpttrack_ = iConfig.getParameter<double> ("MinPtTrack");
17  maxptnearby_ = iConfig.getParameter<double> ("MaxPtNearby");
18  maxetatrack_ = iConfig.getParameter<double> ("MaxEtaTrack");
19  minetatrack_ = iConfig.getParameter<double> ("MinEtaTrack");
20  filterE_ = iConfig.getParameter<bool> ("filterTrackEnergy");
21  minEnergy_ = iConfig.getParameter<double>("MinEnergyTrack");
22  nMaxTrackCandidates_ = iConfig.getParameter<int>("NMaxTrackCandidates");
23  dropMultiL2Event_ = iConfig.getParameter<bool> ("DropMultiL2Event");
24  candToken_ = consumes<reco::IsolatedPixelTrackCandidateCollection>(candTag_);
25  hltGTseedToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(hltGTseedlabel_);
26 }
27 
29 
30 
31 void
35  desc.add<edm::InputTag>("candTag",edm::InputTag("isolPixelTrackProd"));
36  desc.add<edm::InputTag>("L1GTSeedLabel",edm::InputTag("hltL1sIsoTrack"));
37  desc.add<double>("MaxPtNearby",2.0);
38  desc.add<double>("MinEnergyTrack",12.0);
39  desc.add<double>("MinPtTrack",3.5);
40  desc.add<double>("MaxEtaTrack",1.15);
41  desc.add<double>("MinEtaTrack",0.0);
42  desc.add<double>("MinDeltaPtL1Jet",-40000.0);
43  desc.add<bool>("filterTrackEnergy",true);
44  desc.add<int>("NMaxTrackCandidates",10);
45  desc.add<bool>("DropMultiL2Event",false);
46  descriptions.add("hltPixelIsolTrackFilter",desc);
47 }
48 
50 
51  if (saveTags())
52  filterproduct.addCollectionTag(candTag_);
53 
54  // Ref to Candidate object to be recorded in filter object
56 
57  // get hold of filtered candidates
59  iEvent.getByToken(candToken_,recotrackcands);
60 
61  //Filtering
62 
63  //find leading L1 jet:
64  double ptTriggered = -10;
65 
67  iEvent.getByToken(hltGTseedToken_, l1trigobj);
68 
69  std::vector< edm::Ref<l1extra::L1JetParticleCollection> > l1tauobjref;
70  std::vector< edm::Ref<l1extra::L1JetParticleCollection> > l1jetobjref;
71  std::vector< edm::Ref<l1extra::L1JetParticleCollection> > l1forjetobjref;
72 
73  l1trigobj->getObjects(trigger::TriggerL1TauJet, l1tauobjref);
74  l1trigobj->getObjects(trigger::TriggerL1CenJet, l1jetobjref);
75  l1trigobj->getObjects(trigger::TriggerL1ForJet, l1forjetobjref);
76 
77  for (auto & p : l1tauobjref)
78  if (p->pt() > ptTriggered)
79  ptTriggered = p->pt();
80  for (auto & p : l1jetobjref)
81  if (p->pt() > ptTriggered)
82  ptTriggered = p->pt();
83  for (auto & p : l1forjetobjref)
84  if (p->pt() > ptTriggered)
85  ptTriggered = p->pt();
86 
87  int n=0;
88  for (unsigned int i=0; i<recotrackcands->size(); i++) {
90 
91  // cut on deltaPT
92  if (ptTriggered-candref->pt()<minDeltaPtL1Jet_) continue;
93 
94  // select on transverse momentum
95  if (!filterE_&&(candref->maxPtPxl()<maxptnearby_)&&
96  (candref->pt()>minpttrack_)&&fabs(candref->track()->eta())<maxetatrack_&&fabs(candref->track()->eta())>minetatrack_) {
97  filterproduct.addObject(trigger::TriggerTrack, candref);
98  n++;
99  LogDebug("IsoTrk") << "PixelIsolP:Candidate[" << n <<"] pt|eta|phi "
100  << candref->pt() << "|" << candref->eta() << "|"
101  << candref->phi() << "\n";
102  }
103 
104  // select on momentum
105  if (filterE_){
106  if ((candref->maxPtPxl()<maxptnearby_)&&((candref->pt())*cosh(candref->track()->eta())>minEnergy_)&&fabs(candref->track()->eta())<maxetatrack_&&fabs(candref->track()->eta())>minetatrack_) {
107  filterproduct.addObject(trigger::TriggerTrack, candref);
108  n++;
109  LogDebug("IsoTrk") << "PixelIsolE:Candidate[" << n <<"] pt|eta|phi "
110  << candref->pt() << "|" << candref->eta() << "|"
111  << candref->phi() << "\n";
112  }
113  }
114 
115  // stop looping over tracks if max number is reached
116  if(!dropMultiL2Event_ && n>=nMaxTrackCandidates_) break;
117 
118  } // loop over tracks
119 
120 
121  bool accept(n>0);
122 
123  if( dropMultiL2Event_ && n>nMaxTrackCandidates_ ) accept=false;
124 
125  return accept;
126 
127 }
128 
129 
130 // declare this class as a framework plugin
#define LogDebug(id)
T getParameter(std::string const &) const
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
edm::EDGetTokenT< reco::IsolatedPixelTrackCandidateCollection > candToken_
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
~HLTPixelIsolTrackFilter() override
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > hltGTseedToken_
bool saveTags() const
Definition: HLTFilter.h:45
HLTPixelIsolTrackFilter(const edm::ParameterSet &)