CMS 3D CMS Logo

HLTCATopTagFilter.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HLTCATopTagFilter
4 // Class: HLTCATopTagFilter
5 //
13 //
14 // Original Author: dylan rankin
15 // Created: Wed, 17 Jul 2013 22:11:30 GMT
16 // $Id$
17 //
18 //
19 
21 
22 using namespace std;
23 using namespace reco;
24 using namespace edm;
25 
26 //
27 // constructors and destructor
28 //
29 
31  src_ (iConfig.getParameter<edm::InputTag>("src")),
32  pfsrc_ (iConfig.getParameter<edm::InputTag>("pfsrc")),
33  inputToken_ (consumes<reco::BasicJetCollection>(src_)),
34  inputPFToken_ (consumes<reco::PFJetCollection>(pfsrc_))
35 {
36  TopMass_ = iConfig.getParameter<double>("TopMass");
37  minTopMass_ = iConfig.getParameter<double>("minTopMass");
38  maxTopMass_ = iConfig.getParameter<double>("maxTopMass");
39  minMinMass_ = iConfig.getParameter<double>("minMinMass");
40 }
41 
42 
44 
45 
49  desc.add<double>("TopMass",171.);
50  desc.add<double>("maxTopMass",230.);
51  desc.add<double>("minTopMass",140.);
52  desc.add<double>("minMinMass",50.);
53  desc.add<edm::InputTag>("src",edm::InputTag("hltParticleFlow"));
54  desc.add<edm::InputTag>("pfsrc",edm::InputTag("selectedPFJets"));
55  desc.add<int>("triggerType",trigger::TriggerJet);
56  descriptions.add("hltCA8TopTagFilter",desc);
57 }
58 // ------------ method called to for each event ------------
59 
61 {
62 
63 
64  //get basic jets
66  iEvent.getByToken(inputToken_, pBasicJets);
67 
68  //get corresponding pf jets
70  iEvent.getByToken(inputPFToken_, pfJets);
71 
72  //add filter object
73  if(saveTags()){
74  filterobject.addCollectionTag(pfsrc_);
75  }
76 
77  //initialize the properties
79  CATopJetProperties properties;
80 
81  // Now loop over the hard jets and do kinematic cuts
82  auto ihardJet = pBasicJets->begin(),
83  ihardJetEnd = pBasicJets->end();
84  auto ipfJet = pfJets->begin();
85  bool pass = false;
86 
87  for ( ; ihardJet != ihardJetEnd; ++ihardJet, ++ipfJet ) {
88 
89  //if (ihardJet->pt() < 350) continue;
90 
91  // Get properties
92  properties = helper( (reco::Jet&) *ihardJet );
93 
94  if (properties.nSubJets < 3 ||properties.minMass < minMinMass_ || properties.topMass < minTopMass_ || properties.topMass > maxTopMass_) continue;
95  else {
96  // Get a ref to the hard jet
97  reco::PFJetRef ref = reco::PFJetRef(pfJets,distance(pfJets->begin(),ipfJet));
98  //add ref to event
99  filterobject.addObject(trigger::TriggerJet,ref);
100  pass = true;
101  }
102 
103  }// end loop over hard jets
104 
105 
106 
107 
108  return pass;
109 }
110 // ------------ method called once each job just before starting event loop ------------
111 
112 
113 
114 // declare this class as a framework plugin
T getParameter(std::string const &) const
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterobject) const override
Definition: helper.py:1
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
Base class for all types of Jets.
Definition: Jet.h:20
edm::Ref< PFJetCollection > PFJetRef
edm references
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ParameterDescriptionBase * add(U const &iLabel, T const &value)
const edm::EDGetTokenT< reco::BasicJetCollection > inputToken_
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)
HLTCATopTagFilter(const edm::ParameterSet &)
std::vector< PFJet > PFJetCollection
collection of PFJet objects
fixed size matrix
bool saveTags() const
Definition: HLTFilter.h:45
HLT enums.
const edm::EDGetTokenT< reco::PFJetCollection > inputPFToken_
~HLTCATopTagFilter() override
edm::InputTag pfsrc_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::vector< BasicJet > BasicJetCollection
collection of BasicJet objects