CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTCaloJetIDProducer.cc
Go to the documentation of this file.
4 
8 
10  jetsInput_ (iConfig.getParameter<edm::InputTag>("jetsInput")),
11  min_EMF_ (iConfig.getParameter<double>("min_EMF")),
12  max_EMF_ (iConfig.getParameter<double>("max_EMF")),
13  min_N90_ (iConfig.getParameter<int>("min_N90")),
14  min_N90hits_ (iConfig.getParameter<int>("min_N90hits")),
15  jetID_ (iConfig.getParameter<edm::ParameterSet>("JetIDParams"))
16 {
17  // produces< reco::CaloJetCollection > ( "hltCaloJetIDCollection" );
18  produces< reco::CaloJetCollection > ();
19  m_theCaloJetToken = consumes<reco::CaloJetCollection>(jetsInput_);
20 }
21 
23 {
24 
25 }
26 
28 {
29 
30 }
31 
32 void
35  desc.add<edm::InputTag>("jetsInput",edm::InputTag("hltMCJetCorJetIcone5HF07"));
36  desc.add<double>("min_EMF",0.0001);
37  desc.add<double>("max_EMF",999.);
38  desc.add<int>("min_N90",0);
39  desc.add<int>("min_N90hits",2);
41  jetidPSet.add<bool>("useRecHits",true);
42  jetidPSet.add<edm::InputTag>("hbheRecHitsColl",edm::InputTag("hltHbhereco"));
43  jetidPSet.add<edm::InputTag>("hoRecHitsColl",edm::InputTag("hltHoreco"));
44  jetidPSet.add<edm::InputTag>("hfRecHitsColl",edm::InputTag("hltHfreco"));
45  jetidPSet.add<edm::InputTag>("ebRecHitsColl",edm::InputTag("hltEcalRecHitAll", "EcalRecHitsEB"));
46  jetidPSet.add<edm::InputTag>("eeRecHitsColl",edm::InputTag("hltEcalRecHitAll", "EcalRecHitsEE"));
47  desc.add<edm::ParameterSetDescription>("JetIDParams",jetidPSet);
48  descriptions.add("hltCaloJetIDProducer",desc);
49 }
50 
51 
53 {
54 
56  iEvent.getByToken(m_theCaloJetToken, calojets);
57 
58  std::auto_ptr<reco::CaloJetCollection> result (new reco::CaloJetCollection);
59 
60  for (reco::CaloJetCollection::const_iterator calojetc = calojets->begin();
61  calojetc != calojets->end(); ++calojetc) {
62 
63  if (std::abs(calojetc->eta()) >= 2.6) {
64  result->push_back(*calojetc);
65  } else {
66  if (min_N90hits_>0) jetID_.calculate( iEvent, *calojetc );
67  if ((calojetc->emEnergyFraction() >= min_EMF_) && ((min_N90hits_<=0) || (jetID_.n90Hits() >= min_N90hits_)) && (calojetc->n90() >= min_N90_) && (calojetc->emEnergyFraction() <= max_EMF_)) {
68  result->push_back(*calojetc);
69  }
70  }
71  } // calojetc
72 
73  //iEvent.put( result, "hltCaloJetIDCollection");
74  iEvent.put( result);
75 
76 }
HLTCaloJetIDProducer(const edm::ParameterSet &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
virtual void produce(edm::Event &, const edm::EventSetup &)
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
tuple result
Definition: query.py:137
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 add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< reco::CaloJetCollection > m_theCaloJetToken
reco::helper::JetIDHelper jetID_
void calculate(const edm::Event &event, const reco::CaloJet &jet, const int iDbg=0)
Definition: JetIDHelper.cc:86
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects