CMS 3D CMS Logo

HLTJetCollectionsVBFFilter.cc
Go to the documentation of this file.
1 
19 
20 
21 //
22 // constructors and destructor
23 //
24 template <typename T>
26  inputTag_(iConfig.getParameter< edm::InputTag > ("inputTag")),
27  originalTag_(iConfig.getParameter< edm::InputTag > ("originalTag")),
28  softJetPt_(iConfig.getParameter<double> ("SoftJetPt")),
29  hardJetPt_(iConfig.getParameter<double> ("HardJetPt")),
30  minDeltaEta_(iConfig.getParameter<double> ("MinDeltaEta")),
31  thirdJetPt_(iConfig.getParameter<double> ("ThirdJetPt")),
32  maxAbsJetEta_(iConfig.getParameter<double> ("MaxAbsJetEta")),
33  maxAbsThirdJetEta_(iConfig.getParameter<double> ("MaxAbsThirdJetEta")),
34  minNJets_(iConfig.getParameter<unsigned int> ("MinNJets")),
35  triggerType_(iConfig.getParameter<int> ("TriggerType"))
36 {
37  typedef std::vector<edm::RefVector<std::vector<T>,T,edm::refhelper::FindUsingAdvance<std::vector<T>,T> > > TCollectionVector;
38  m_theJetToken = consumes<TCollectionVector>(inputTag_);
39 }
40 
41 template <typename T>
43 
44 template <typename T>
45 void
49  desc.add<edm::InputTag>("inputTag",edm::InputTag("hltIterativeCone5CaloJets"));
50  desc.add<edm::InputTag>("originalTag",edm::InputTag("hltIterativeCone5CaloJets"));
51  desc.add<double>("SoftJetPt",25.0);
52  desc.add<double>("HardJetPt",35.0);
53  desc.add<double>("MinDeltaEta",3.0);
54  desc.add<double>("ThirdJetPt",20.0);
55  desc.add<double>("MaxAbsJetEta",9999.);
56  desc.add<double>("MaxAbsThirdJetEta",2.6);
57  desc.add<unsigned int>("MinNJets",2);
58  desc.add<int>("TriggerType",trigger::TriggerJet);
60 }
61 
62 // ------------ method called to produce the data ------------
63 template <typename T>
64 bool
66 {
67  using namespace std;
68  using namespace edm;
69  using namespace reco;
70  using namespace trigger;
71 
72  typedef vector<T> TCollection;
73  typedef Ref<TCollection> TRef;
74  typedef edm::RefVector<TCollection> TRefVector;
75  typedef std::vector<edm::RefVector<std::vector<T>,T,edm::refhelper::FindUsingAdvance<std::vector<T>,T> > > TCollectionVector;
76 
77  // The filter object
78  if (saveTags()) filterproduct.addCollectionTag(originalTag_);
79 
80  Handle<TCollectionVector> theJetCollectionsHandle;
81  iEvent.getByToken(m_theJetToken, theJetCollectionsHandle);
82  const TCollectionVector & theJetCollections = *theJetCollectionsHandle;
83  // filter decision
84  bool accept(false);
85  std::vector < TRef > goodJetRefs;
86 
87  for(unsigned int collection = 0; collection < theJetCollections.size(); ++ collection) {
88 
89  const TRefVector & refVector = theJetCollections[collection];
90  if(refVector.size() < minNJets_) continue;
91 
92  // VBF decision
93  bool thereAreVBFJets(false);
94  // 3rd Jet check decision
95  bool goodThirdJet(false);
96  if ( minNJets_ < 3 ) goodThirdJet = true;
97 
98  //empty the good jets collection
99  goodJetRefs.clear();
100 
101  TRef refOne;
102  TRef refTwo;
103  typename TRefVector::const_iterator jetOne ( refVector.begin() );
104  int firstJetIndex=100, secondJetIndex=100, thirdJetIndex=100;
105 
106  // Cycle to look for VBF jets
107  for (; jetOne != refVector.end(); jetOne++) {
108  TRef jetOneRef(*jetOne);
109 
110  if ( thereAreVBFJets ) break;
111  if ( jetOneRef->pt() < hardJetPt_ ) break;
112  if ( std::abs(jetOneRef->eta()) > maxAbsJetEta_ ) continue;
113 
114  typename TRefVector::const_iterator jetTwo = jetOne + 1;
115  secondJetIndex = firstJetIndex;
116  for (; jetTwo != refVector.end(); jetTwo++) {
117  TRef jetTwoRef(*jetTwo);
118 
119  if ( jetTwoRef->pt() < softJetPt_ ) break;
120  if ( std::abs(jetTwoRef->eta()) > maxAbsJetEta_ ) continue;
121 
122  if ( std::abs(jetTwoRef->eta() - jetOneRef->eta()) < minDeltaEta_ ) continue;
123 
124  thereAreVBFJets = true;
125  refOne = *jetOne;
126  goodJetRefs.push_back(refOne);
127  refTwo = *jetTwo;
128  goodJetRefs.push_back(refTwo);
129 
130  firstJetIndex = (int) (jetOne - refVector.begin());
131  secondJetIndex= (int) (jetTwo - refVector.begin());
132 
133  break;
134 
135  }
136  }// Close looop on VBF
137 
138  // Look for a third jet, if you've found the previous 2
139  if ( minNJets_ > 2 && thereAreVBFJets ) {
140  TRef refThree;
141  typename TRefVector::const_iterator jetThree ( refVector.begin() );
142  for (; jetThree != refVector.end(); jetThree++) {
143  thirdJetIndex = (int) (jetThree - refVector.begin());
144 
145  TRef jetThreeRef(*jetThree);
146 
147  if ( thirdJetIndex == firstJetIndex || thirdJetIndex == secondJetIndex ) continue;
148 
149  if (jetThreeRef->pt() >= thirdJetPt_ && std::abs(jetThreeRef->eta()) <= maxAbsThirdJetEta_) {
150  goodThirdJet = true;
151  refThree = *jetThree;
152  goodJetRefs.push_back(refThree);
153  break;
154  }
155  }
156  }
157 
158  if(thereAreVBFJets && goodThirdJet){
159  accept = true;
160  break;
161  }
162 
163  }
164 
165  //fill the filter object
166  for (unsigned int refIndex = 0; refIndex < goodJetRefs.size(); ++refIndex) {
167  filterproduct.addObject(triggerType_, goodJetRefs.at(refIndex));
168  }
169 
170  return accept;
171 }
edm::EDGetTokenT< std::vector< edm::RefVector< std::vector< T >, T, edm::refhelper::FindUsingAdvance< std::vector< T >, T > > > > m_theJetToken
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
~HLTJetCollectionsVBFFilter() override
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
int iEvent
Definition: GenABIO.cc:224
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
HLTJetCollectionsVBFFilter(const edm::ParameterSet &)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
bool saveTags() const
Definition: HLTFilter.h:45
HLT enums.
long double T