CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTJetVBFFilter.cc
Go to the documentation of this file.
1 
11 
15 
23 
24 #include<typeinfo>
25 
26 //
27 // constructors and destructor
28 //
29 template<typename T>
31  inputTag_ (iConfig.template getParameter< edm::InputTag > ("inputTag")),
32  minPtLow_ (iConfig.template getParameter<double> ("minPtLow")),
33  minPtHigh_ (iConfig.template getParameter<double> ("minPtHigh")),
34  etaOpposite_ (iConfig.template getParameter<bool> ("etaOpposite")),
35  minDeltaEta_ (iConfig.template getParameter<double> ("minDeltaEta")),
36  minInvMass_ (iConfig.template getParameter<double> ("minInvMass")),
37  maxEta_ (iConfig.template getParameter<double> ("maxEta")),
38  leadingJetOnly_ (iConfig.template getParameter<bool> ("leadingJetOnly")),
39  triggerType_ (iConfig.template getParameter<int> ("triggerType"))
40 {
41  LogDebug("") << "HLTJetVBFFilter: Input/minPtLow_/minPtHigh_/triggerType : "
42  << inputTag_.encode() << " "
43  << minPtLow_ << " "
44  << minPtHigh_ << " "
45  << triggerType_;
46 }
47 
48 template<typename T>
50 
51 template<typename T>
52 void
55  makeHLTFilterDescription(desc);
56  desc.add<edm::InputTag>("inputTag",edm::InputTag("hltAntiKT5ConvPFJets"));
57  desc.add<double>("minPtLow",40.);
58  desc.add<double>("minPtHigh",40.);
59  desc.add<bool>("etaOpposite",false);
60  desc.add<double>("minDeltaEta",4.0);
61  desc.add<double>("minInvMass",1000.);
62  desc.add<double>("maxEta",5.0);
63  desc.add<bool>("leadingJetOnly",false);
64  desc.add<int>("triggerType",trigger::TriggerJet);
65  descriptions.add(std::string("hlt")+std::string(typeid(HLTJetVBFFilter<T>).name()),desc);
66 }
67 
68 //
69 // ------------ method called to produce the data ------------
70 //
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 
83  // The filter object
84  if (saveTags()) filterproduct.addCollectionTag(inputTag_);
85 
86  // get hold of collection of objects
87  Handle<TCollection> objects;
88  iEvent.getByLabel (inputTag_,objects);
89 
90  // look at all candidates, check cuts and add to filter object
91  int n(0);
92 
93  // events with two or more jets
94  if(objects->size() > 1){
95 
96  double ejet1 = 0.;
97  double pxjet1 = 0.;
98  double pyjet1 = 0.;
99  double pzjet1 = 0.;
100  double ptjet1 = 0.;
101  double etajet1 = 0.;
102 
103  double ejet2 = 0.;
104  double pxjet2 = 0.;
105  double pyjet2 = 0.;
106  double pzjet2 = 0.;
107  double ptjet2 = 0.;
108  double etajet2 = 0.;
109 
110  // loop on all jets
111  int countJet1(0);
112  int countJet2(0);
113  typename TCollection::const_iterator jet1 ( objects->begin() );
114  for (; jet1!=objects->end(); jet1++) {
115  countJet1++;
116  if( leadingJetOnly_==true && countJet1>2 ) break;
117  //
118  if( jet1->pt() < minPtHigh_ ) break; //No need to go to the next jet (lower PT)
119  if( std::abs(jet1->eta()) > maxEta_ ) continue;
120  //
121  countJet2 = countJet1-1;
122  typename TCollection::const_iterator jet2 ( jet1+1 );
123  for (; jet2!=objects->end(); jet2++) {
124  countJet2++;
125  if( leadingJetOnly_==true && countJet2>2 ) break;
126  //
127  if( jet2->pt() < minPtLow_ ) break; //No need to go to the next jet (lower PT)
128  if( std::abs(jet2->eta()) > maxEta_ ) continue;
129  //
130  ejet1 = jet1->energy();
131  pxjet1 = jet1->px();
132  pyjet1 = jet1->py();
133  pzjet1 = jet1->pz();
134  ptjet1 = jet1->pt();
135  etajet1 = jet1->eta();
136 
137  ejet2 = jet2->energy();
138  pxjet2 = jet2->px();
139  pyjet2 = jet2->py();
140  pzjet2 = jet2->pz();
141  ptjet2 = jet2->pt();
142  etajet2 = jet2->eta();
143  //
144  float deltaetajet = etajet1 - etajet2;
145  float invmassjet = sqrt( (ejet1 + ejet2) * (ejet1 + ejet2) -
146  (pxjet1 + pxjet2) * (pxjet1 + pxjet2) -
147  (pyjet1 + pyjet2) * (pyjet1 + pyjet2) -
148  (pzjet1 + pzjet2) * (pzjet1 + pzjet2) );
149 
150  // VBF cuts
151  if ( (ptjet1 > minPtHigh_) &&
152  (ptjet2 > minPtLow_) &&
153  ( (etaOpposite_ == true && etajet1*etajet2 < 0) || (etaOpposite_ == false) ) &&
154  (std::abs(deltaetajet) > minDeltaEta_) &&
155  (std::abs(invmassjet) > minInvMass_) ){
156  ++n;
157  TRef ref1 = TRef(objects,distance(objects->begin(),jet1));
158  TRef ref2 = TRef(objects,distance(objects->begin(),jet2));
159  filterproduct.addObject(triggerType_,ref1);
160  filterproduct.addObject(triggerType_,ref2);
161  }// VBF cuts
162  //if(n>=1) break; //Store all possible pairs
163  }
164  //if(n>=1) break; //Store all possible pairs
165  }// loop on all jets
166  }// events with two or more jets
167 
168  // filter decision
169  bool accept(n>=1);
170 
171  return accept;
172 }
#define LogDebug(id)
#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
std::string encode() const
Definition: InputTag.cc:72
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
edm::InputTag inputTag_
T sqrt(T t)
Definition: SSEVec.h:46
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)
HLTJetVBFFilter(const edm::ParameterSet &)
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
def template
Definition: svgfig.py:520