CMS 3D CMS Logo

HLTForwardBackwardJetsFilter.cc
Go to the documentation of this file.
1 
8 
12 
21 
22 //
23 // constructors and destructor
24 //
25 template<typename T>
27  HLTFilter(iConfig),
28  inputTag_ (iConfig.template getParameter< edm::InputTag > ("inputTag")),
29  minPt_ (iConfig.template getParameter<double> ("minPt")),
30  minEta_ (iConfig.template getParameter<double> ("minEta")),
31  maxEta_ (iConfig.template getParameter<double> ("maxEta")),
32  nNeg_ (iConfig.template getParameter<unsigned int>("nNeg")),
33  nPos_ (iConfig.template getParameter<unsigned int>("nPos")),
34  nTot_ (iConfig.template getParameter<unsigned int>("nTot")),
35  triggerType_ (iConfig.template getParameter<int> ("triggerType"))
36 {
37  m_theJetToken = consumes<std::vector<T>>(inputTag_);
38  LogDebug("") << "HLTForwardBackwardJetsFilter: Input/minPt/minEta/maxEta/triggerType : "
39  << inputTag_.encode() << " "
40  << minPt_ << " "
41  << minEta_ << " "
42  << maxEta_ << " "
43  << nNeg_ << " "
44  << nPos_ << " "
45  << nTot_ << " "
46  << triggerType_;
47 }
48 
49 template<typename T>
51 
52 template<typename T>
53 void
57  desc.add<edm::InputTag>("inputTag",edm::InputTag("hltIterativeCone5CaloJetsRegional"));
58  desc.add<double>("minPt",15.0);
59  desc.add<double>("minEta",3.0);
60  desc.add<double>("maxEta",5.1);
61  desc.add<unsigned int>("nNeg",1);
62  desc.add<unsigned int>("nPos",1);
63  desc.add<unsigned int>("nTot",0);
64  desc.add<int>("triggerType",trigger::TriggerJet);
66 }
67 
68 // ------------ method called to produce the data ------------
69 template<typename T>
70 bool
72 {
73  using namespace std;
74  using namespace edm;
75  using namespace reco;
76  using namespace trigger;
77 
78  typedef vector<T> TCollection;
79  typedef Ref<TCollection> TRef;
80 
81  // The filter object
82  if (saveTags()) filterproduct.addCollectionTag(inputTag_);
83 
84  // get hold of collection of objects
86  iEvent.getByToken(m_theJetToken,objects);
87 
88  // look at all candidates, check cuts and add to filter object
89  unsigned int nPosJets(0);
90  unsigned int nNegJets(0);
91 
92  typename TCollection::const_iterator jet;
93  // look for jets satifying pt and eta cuts; first on the plus side, then the minus side
94 
95  for (jet=objects->begin(); jet!=objects->end(); jet++) {
96  double ptjet = jet->pt();
97  double etajet = jet->eta();
98  if( ptjet >= minPt_ ){
99  if (( minEta_<= etajet) && (etajet <= maxEta_) ){
100  nPosJets++;
101  TRef ref = TRef(objects,distance(objects->begin(),jet));
102  filterproduct.addObject(triggerType_,ref);
103  }
104  if ((-maxEta_<= etajet) && (etajet <=-minEta_) ){
105  nNegJets++;
106  TRef ref = TRef(objects,distance(objects->begin(),jet));
107  filterproduct.addObject(triggerType_,ref);
108  }
109  }
110  }
111 
112  // filter decision
113  const bool accept(
114  ( nNegJets >= nNeg_ ) &&
115  ( nPosJets >= nPos_ ) &&
116  ((nNegJets+nPosJets) >= nTot_ )
117  );
118 
119  return accept;
120 }
#define LogDebug(id)
~HLTForwardBackwardJetsFilter() override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
std::string defaultModuleLabel()
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
std::string encode() const
Definition: InputTag.cc:159
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
int iEvent
Definition: GenABIO.cc:224
def template(fileName, svg, replaceme="REPLACEME")
Definition: svgfig.py:521
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
fixed size matrix
bool saveTags() const
Definition: HLTFilter.h:45
HLT enums.
HLTForwardBackwardJetsFilter(const edm::ParameterSet &)
edm::EDGetTokenT< std::vector< T > > m_theJetToken