CMS 3D CMS Logo

HLTDiJetAveEtaFilter.cc
Go to the documentation of this file.
1 
20 
21 
22 //
23 // constructors and destructor
24 //
25 template<typename T>
27  inputJetTag_ (iConfig.template getParameter< edm::InputTag > ("inputJetTag")),
28  minPtJet_ (iConfig.template getParameter<double> ("minPtJet")),
29  minPtAve_ (iConfig.template getParameter<double> ("minPtAve")),
30  //minPtJet3_ (iConfig.template getParameter<double> ("minPtJet3")),
31  minDphi_ (iConfig.template getParameter<double> ("minDphi")),
32  tagEtaMin_ (iConfig.template getParameter<double> ("minTagEta")),
33  tagEtaMax_ (iConfig.template getParameter<double> ("maxTagEta")),
34  probeEtaMin_ (iConfig.template getParameter<double> ("minProbeEta")),
35  probeEtaMax_ (iConfig.template getParameter<double> ("maxProbeEta")),
36  triggerType_ (iConfig.template getParameter<int> ("triggerType"))
37 {
38  m_theJetToken = consumes<std::vector<T>>(inputJetTag_);
39  LogDebug("") << "HLTDiJetAveEtaFilter: Input/minPtAve/minDphi/triggerType : "
40  << inputJetTag_.encode() << " "
41  << minPtAve_ << " "
42  //<< minPtJet3_ << " "
43  << minDphi_ << " "
44  << triggerType_;
45 }
46 
47 template<typename T>
49 
50 template<typename T>
51 void
55  desc.add<edm::InputTag>("inputJetTag",edm::InputTag("hltIterativeCone5CaloJets"));
56  desc.add<double>("minPtAve",100.0);
57  desc.add<double>("minPtJet",50.0);
58  //desc.add<double>("minPtJet3",99999.0);
59  desc.add<double>("minDphi",-1.0);
60  desc.add<double>("minTagEta", -1.);
61  desc.add<double>("maxTagEta", 1.4);
62  desc.add<double>("minProbeEta", 2.7);
63  desc.add<double>("maxProbeEta", 5.5);
64  desc.add<int>("triggerType",trigger::TriggerJet);
65  descriptions.add(defaultModuleLabel<HLTDiJetAveEtaFilter<T>>(), desc);
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(inputJetTag_);
83 
84  // get hold of collection of objects
86  iEvent.getByToken (m_theJetToken,objects);
87 
88  int n(0);
89 
90  if(objects->size() > 1){ // events with two or more jets
91  typename TCollection::const_iterator iProbe ( objects->begin() );
92  typename TCollection::const_iterator iEnd ( objects->end() );
93  for (; iProbe!=iEnd; ++iProbe) {
94  if (iProbe->pt() < minPtJet_) continue;
95 
96  // for easier trigger efficiency evaluation save all probe/tag
97  // objects passing the minPT/eta criteria (outer loop)
98  float eta = std::abs(iProbe->eta());
99  bool isGood = false; // probe or tag
100  bool isProbe = false;
101  if ( eta > probeEtaMin_ && eta < probeEtaMax_ ){
102  isGood = true;
103  isProbe = true;
104  }
105  if ( eta > tagEtaMin_ && eta < tagEtaMax_ ){
106  isGood = true;
107  }
108  if (isGood){
109  filterproduct.addObject(triggerType_, TRef(objects,distance(objects->begin(),iProbe)));
110  }
111 
112  if (!isProbe) continue;
113 
114  typename TCollection::const_iterator iTag ( objects->begin() );
115  for (;iTag != iEnd; ++iTag){
116  if (iTag==iProbe) continue;
117  if (iTag->pt() < minPtJet_) continue;
118  float eta2 = std::abs(iTag->eta());
119  if ( eta2 < tagEtaMin_ || eta2 > tagEtaMax_ ) continue;
120  double dphi = std::abs(deltaPhi(iProbe->phi(),iTag->phi() ));
121  if (dphi<minDphi_) {
122  continue;
123  }
124 
125  double ptAve = (iProbe->pt() + iTag->pt())/2;
126  if (ptAve<minPtAve_ ) {
127  continue;
128  }
129  ++n;
130  }
131  }
132  } // events with two or more jets
133  // filter decision
134  bool accept(n>=1);
135  return accept;
136 }
#define LogDebug(id)
edm::EDGetTokenT< std::vector< 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
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
HLTDiJetAveEtaFilter(const edm::ParameterSet &)
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
~HLTDiJetAveEtaFilter() override
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
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.