CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTFatJetMassFilter.cc
Go to the documentation of this file.
1 
9 
14 
17 
25 
26 #include <vector>
27 
28 #include <typeinfo>
29 
30 
31 //
32 // constructors and destructor
33 //
34 template<typename jetType>
36  HLTFilter(iConfig),
37  inputJetTag_ (iConfig.template getParameter< edm::InputTag > ("inputJetTag")),
38  minMass_ (iConfig.template getParameter<double> ("minMass")),
39  fatJetDeltaR_ (iConfig.template getParameter<double> ("fatJetDeltaR")),
40  maxDeltaEta_ (iConfig.template getParameter<double> ("maxDeltaEta")),
41  maxJetEta_ (iConfig.template getParameter<double> ("maxJetEta")),
42  minJetPt_ (iConfig.template getParameter<double> ("minJetPt")),
43  triggerType_ (iConfig.template getParameter<int> ("triggerType"))
44 {
45  LogDebug("") << "HLTFatJetMassFilter: Input/minMass/fatJetDeltaR/maxDeltaEta/maxJetEta/minJetPt/triggerType : "
46  << inputJetTag_.encode() << " "
47  << minMass_ << " "
48  << fatJetDeltaR_ << " "
49  << maxDeltaEta_ << " "
50  << maxJetEta_ << " "
51  << minJetPt_ << " "
52  << triggerType_;
53 }
54 
55 template<typename jetType>
57 
58 template<typename jetType>
59 void
62  makeHLTFilterDescription(desc);
63  desc.add<edm::InputTag>("inputJetTag",edm::InputTag("hltCollection"));
64  desc.add<double>("minMass",0.0);
65  desc.add<double>("fatJetDeltaR",1.1);
66  desc.add<double>("maxDeltaEta",10.0);
67  desc.add<double>("maxJetEta",3.0);
68  desc.add<double>("minJetPt",30.0);
69  desc.add<int>("triggerType",trigger::TriggerJet);
70  descriptions.add(std::string("hlt")+std::string(typeid(HLTFatJetMassFilter<jetType>).name()),desc);
71 }
72 
73 // ------------ method called to produce the data ------------
74 template<typename jetType>
75 bool
77 {
78  using namespace std;
79  using namespace edm;
80  using namespace reco;
81 
82  typedef vector<jetType> JetCollection;
83  typedef Ref<JetCollection> JetRef;
84 
85  // The filter object
86  if (saveTags()) filterproduct.addCollectionTag(inputJetTag_);
87 
88  // All jets
89  Handle<JetCollection> objects;
90  iEvent.getByLabel (inputJetTag_,objects);
91 
92  // Selected jets
93  CaloJetCollection recojets;
94  typename JetCollection::const_iterator i ( objects->begin() );
95  for(;i != objects->end(); i++){
96  if(std::abs(i->eta()) < maxJetEta_ && i->pt() >= minJetPt_){
97  reco::CaloJet jet(i->p4(), i->vertex(), reco::CaloJet::Specific());
98  recojets.push_back(jet);
99  }
100  }
101 
102  // events with at least two jets
103  if(recojets.size() < 2) return false;
104 
105  math::PtEtaPhiMLorentzVector j1(0.1, 0., 0., 0.);
106  math::PtEtaPhiMLorentzVector j2(0.1, 0., 0., 0.);
107  double jetPt1 = 0.;
108  double jetPt2 = 0.;
109 
110  // look for the two highest-pT jet
111  CaloJetCollection::const_iterator recojet ( recojets.begin() );
112  for (; recojet != recojets.end() ; recojet++) {
113  if(recojet->pt() > jetPt1) {
114  // downgrade the 1st jet to 2nd jet
115  j2 = j1;
116  jetPt2 = j2.pt();
117  // promote this jet to 1st jet
118  j1 = recojet->p4();
119  jetPt1 = recojet->pt();
120  } else if(recojet->pt() > jetPt2) {
121  // promote this jet to 2nd jet
122  j2 = recojet->p4();
123  jetPt2 = recojet->pt();
124  }
125  }
126 
127  // apply DeltaEta cut
128  double DeltaEta = std::abs(j1.eta() - j2.eta());
129  if(DeltaEta > maxDeltaEta_) return false;
130 
133 
134  // apply radiation recovery
135  for ( recojet = recojets.begin() ; recojet != recojets.end() ; recojet++) {
136  double DeltaR1 = sqrt(pow(recojet->phi()-j1.phi(), 2.)+pow(recojet->eta()-j1.eta(),2.));
137  double DeltaR2 = sqrt(pow(recojet->phi()-j2.phi(), 2.)+pow(recojet->eta()-j2.eta(),2.));
138  if(DeltaR1 < DeltaR2 && DeltaR1 < fatJetDeltaR_) {
139  fj1 += recojet->p4();
140  } else if(DeltaR2 < fatJetDeltaR_) {
141  fj2 += recojet->p4();
142  }
143  }
144 
145  // Apply mass cut
146  fj1 += fj2;
147  if(fj1.mass() < minMass_) return false;
148 
149  return true;
150 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
Jets made from CaloTowers.
Definition: CaloJet.h:30
std::vector< Jet > JetCollection
Definition: Jet.h:50
#define abs(x)
Definition: mlp_lapack.h:159
edm::Ref< JetCollection > JetRef
Definition: Jet.h:52
virtual bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct)
std::string encode() const
Definition: InputTag.cc:72
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:26
int iEvent
Definition: GenABIO.cc:243
T sqrt(T t)
Definition: SSEVec.h:46
edm::InputTag inputJetTag_
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)
HLTFatJetMassFilter(const edm::ParameterSet &)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
def template
Definition: svgfig.py:520
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects