CMS 3D CMS Logo

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