CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTDiJetAveFilter.cc
Go to the documentation of this file.
1 
20 
21 //
22 // constructors and destructor
23 //
24 template<typename T>
26  inputJetTag_ (iConfig.template getParameter< edm::InputTag > ("inputJetTag")),
27  minPtAve_ (iConfig.template getParameter<double> ("minPtAve")),
28  minPtJet3_ (iConfig.template getParameter<double> ("minPtJet3")),
29  minDphi_ (iConfig.template getParameter<double> ("minDphi")),
30  triggerType_ (iConfig.template getParameter<int> ("triggerType"))
31 {
32  m_theJetToken = consumes<std::vector<T>>(inputJetTag_);
33  LogDebug("") << "HLTDiJetAveFilter: Input/minPtAve/minPtJet3/minDphi/triggerType : "
34  << inputJetTag_.encode() << " "
35  << minPtAve_ << " "
36  << minPtJet3_ << " "
37  << minDphi_ << " "
38  << triggerType_;
39 }
40 
41 template<typename T>
43 
44 template<typename T>
45 void
48  makeHLTFilterDescription(desc);
49  desc.add<edm::InputTag>("inputJetTag",edm::InputTag("hltIterativeCone5CaloJets"));
50  desc.add<double>("minPtAve",100.0);
51  desc.add<double>("minPtJet3",99999.0);
52  desc.add<double>("minDphi",-1.0);
53  desc.add<int>("triggerType",trigger::TriggerJet);
54  descriptions.add(defaultModuleLabel<HLTDiJetAveFilter<T>>(), desc);
55 }
56 
57 // ------------ method called to produce the data ------------
58 template<typename T>
59 bool
61 {
62  using namespace std;
63  using namespace edm;
64  using namespace reco;
65  using namespace trigger;
66 
67  typedef vector<T> TCollection;
68  typedef Ref<TCollection> TRef;
69 
70  // The filter object
71  if (saveTags()) filterproduct.addCollectionTag(inputJetTag_);
72 
73  // get hold of collection of objects
74  Handle<TCollection> objects;
75  iEvent.getByToken (m_theJetToken,objects);
76 
77  // look at all candidates, check cuts and add to filter object
78  int n(0);
79 
80  if(objects->size() > 1){
81  // events with two or more jets
82 
83  double ptjet1=0., ptjet2=0.,ptjet3=0.;
84  double phijet1=0.,phijet2=0;
85  int countjets =0;
86 
87  int nmax=1;
88  if (objects->size() > 2) nmax=2;
89 
90  TRef JetRef1,JetRef2;
91 
92  typename TCollection::const_iterator i ( objects->begin() );
93  for (; i<=(objects->begin()+nmax); i++) {
94  if(countjets==0) {
95  ptjet1 = i->pt();
96  phijet1 = i->phi();
97  JetRef1 = TRef(objects,distance(objects->begin(),i));
98  }
99  if(countjets==1) {
100  ptjet2 = i->pt();
101  phijet2 = i->phi();
102  JetRef2 = TRef(objects,distance(objects->begin(),i));
103  }
104  if(countjets==2) {
105  ptjet3 = i->pt();
106  }
107  ++countjets;
108  }
109 
110  double PtAve=(ptjet1 + ptjet2) / 2.;
111  double Dphi = std::abs(deltaPhi(phijet1,phijet2));
112 
113  if( PtAve>minPtAve_ && ptjet3<minPtJet3_ && Dphi>minDphi_){
114  filterproduct.addObject(triggerType_,JetRef1);
115  filterproduct.addObject(triggerType_,JetRef2);
116  ++n;
117  }
118 
119  } // events with two or more jets
120 
121 
122 
123  // filter decision
124  bool accept(n>=1);
125 
126  return accept;
127 }
#define LogDebug(id)
virtual bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
std::string defaultModuleLabel()
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:24
std::string encode() const
Definition: InputTag.cc:164
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:230
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
edm::EDGetTokenT< std::vector< T > > m_theJetToken
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::InputTag inputJetTag_
HLTDiJetAveFilter(const edm::ParameterSet &)
def template
Definition: svgfig.py:520