CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
JetIDProducer Class Reference

#include <RecoJets/JetProducers/plugins/JetIDProducer.cc>

Inheritance diagram for JetIDProducer:
edm::stream::EDProducer<>

Public Member Functions

 JetIDProducer (const edm::ParameterSet &)
 
 ~JetIDProducer () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Private Member Functions

void produce (edm::Event &, const edm::EventSetup &) override
 

Private Attributes

reco::helper::JetIDHelper helper_
 
edm::EDGetTokenT< edm::View< reco::CaloJet > > input_jet_token_
 
reco::helper::JetMuonHitsIDHelper muHelper_
 
edm::InputTag src_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Description: Produces a value map of jet—> jet Id

Implementation: There are two modes: AOD only, in which case only a subset of the info is written, and RECO, when all the info is written. The AOD-only case will be suitable for the "very loose" jet ID, whereas the RECO case will be globally suitable.

Definition at line 44 of file JetIDProducer.h.

Constructor & Destructor Documentation

◆ JetIDProducer()

JetIDProducer::JetIDProducer ( const edm::ParameterSet iConfig)
explicit

Definition at line 17 of file JetIDProducer.cc.

References input_jet_token_, and src_.

18  : src_(iConfig.getParameter<edm::InputTag>("src")),
19  helper_(iConfig, consumesCollector()),
20  muHelper_(iConfig, consumesCollector()) {
21  produces<reco::JetIDValueMap>();
22 
23  input_jet_token_ = consumes<edm::View<reco::CaloJet> >(src_);
24 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
reco::helper::JetIDHelper helper_
Definition: JetIDProducer.h:54
edm::InputTag src_
Definition: JetIDProducer.h:53
reco::helper::JetMuonHitsIDHelper muHelper_
Definition: JetIDProducer.h:55
edm::EDGetTokenT< edm::View< reco::CaloJet > > input_jet_token_
Definition: JetIDProducer.h:57

◆ ~JetIDProducer()

JetIDProducer::~JetIDProducer ( )
override

Definition at line 26 of file JetIDProducer.cc.

26 {}

Member Function Documentation

◆ produce()

void JetIDProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 33 of file JetIDProducer.cc.

References reco::helper::JetIDHelper::approximatefHPD(), reco::helper::JetIDHelper::approximatefRBX(), reco::helper::JetMuonHitsIDHelper::calculate(), reco::helper::JetIDHelper::calculate(), reco::helper::JetIDHelper::fEB(), reco::helper::JetIDHelper::fEE(), reco::helper::JetIDHelper::fHB(), reco::helper::JetIDHelper::fHE(), reco::helper::JetIDHelper::fHFOOT(), reco::helper::JetIDHelper::fHO(), reco::helper::JetIDHelper::fHPD(), trigObjTnPSource_cfi::filler, reco::helper::JetIDHelper::fLong(), reco::helper::JetIDHelper::fLSbad(), reco::helper::JetIDHelper::fRBX(), reco::helper::JetIDHelper::fShort(), reco::helper::JetIDHelper::fSubDetector1(), reco::helper::JetIDHelper::fSubDetector2(), reco::helper::JetIDHelper::fSubDetector3(), reco::helper::JetIDHelper::fSubDetector4(), helper_, reco::helper::JetIDHelper::hitsInN90(), iEvent, input_jet_token_, eostools::move(), muHelper_, reco::helper::JetIDHelper::n90Hits(), reco::helper::JetIDHelper::nECALTowers(), reco::helper::JetIDHelper::nHCALTowers(), BTaggingMonitoring_cff::njets, reco::helper::JetMuonHitsIDHelper::numberOfHits2RPC(), reco::helper::JetMuonHitsIDHelper::numberOfHits3RPC(), reco::helper::JetMuonHitsIDHelper::numberOfHitsRPC(), and reco::helper::JetIDHelper::restrictedEMF().

33  {
34  // get the input jets
36  iEvent.getByToken(input_jet_token_, h_jets);
37 
38  // allocate the jet--->jetid value map
39  auto jetIdValueMap = std::make_unique<reco::JetIDValueMap>();
40  // instantiate the filler with the map
41  reco::JetIDValueMap::Filler filler(*jetIdValueMap);
42 
43  // allocate the vector of ids
44  size_t njets = h_jets->size();
45  std::vector<reco::JetID> ids(njets);
46 
47  // loop over the jets
48  for (edm::View<reco::CaloJet>::const_iterator jetsBegin = h_jets->begin(), jetsEnd = h_jets->end(), ijet = jetsBegin;
49  ijet != jetsEnd;
50  ++ijet) {
51  // get the id from each jet
52  helper_.calculate(iEvent, iSetup, *ijet);
53 
54  muHelper_.calculate(iEvent, iSetup, *ijet);
55 
56  ids[ijet - jetsBegin].fHPD = helper_.fHPD();
57  ids[ijet - jetsBegin].fRBX = helper_.fRBX();
58  ids[ijet - jetsBegin].n90Hits = helper_.n90Hits();
59  ids[ijet - jetsBegin].fSubDetector1 = helper_.fSubDetector1();
60  ids[ijet - jetsBegin].fSubDetector2 = helper_.fSubDetector2();
61  ids[ijet - jetsBegin].fSubDetector3 = helper_.fSubDetector3();
62  ids[ijet - jetsBegin].fSubDetector4 = helper_.fSubDetector4();
63  ids[ijet - jetsBegin].restrictedEMF = helper_.restrictedEMF();
64  ids[ijet - jetsBegin].nHCALTowers = helper_.nHCALTowers();
65  ids[ijet - jetsBegin].nECALTowers = helper_.nECALTowers();
66  ids[ijet - jetsBegin].approximatefHPD = helper_.approximatefHPD();
67  ids[ijet - jetsBegin].approximatefRBX = helper_.approximatefRBX();
68  ids[ijet - jetsBegin].hitsInN90 = helper_.hitsInN90();
69 
70  ids[ijet - jetsBegin].numberOfHits2RPC = muHelper_.numberOfHits2RPC();
71  ids[ijet - jetsBegin].numberOfHits3RPC = muHelper_.numberOfHits3RPC();
72  ids[ijet - jetsBegin].numberOfHitsRPC = muHelper_.numberOfHitsRPC();
73 
74  ids[ijet - jetsBegin].fEB = helper_.fEB();
75  ids[ijet - jetsBegin].fEE = helper_.fEE();
76  ids[ijet - jetsBegin].fHB = helper_.fHB();
77  ids[ijet - jetsBegin].fHE = helper_.fHE();
78  ids[ijet - jetsBegin].fHO = helper_.fHO();
79  ids[ijet - jetsBegin].fLong = helper_.fLong();
80  ids[ijet - jetsBegin].fShort = helper_.fShort();
81  ids[ijet - jetsBegin].fLS = helper_.fLSbad();
82  ids[ijet - jetsBegin].fHFOOT = helper_.fHFOOT();
83  }
84 
85  // set up the map
86  filler.insert(h_jets, ids.begin(), ids.end());
87 
88  // fill the vals
89  filler.fill();
90 
91  // write map to the event
92  iEvent.put(std::move(jetIdValueMap));
93 }
double fSubDetector1() const
Definition: JetIDHelper.h:49
double approximatefRBX() const
Definition: JetIDHelper.h:68
double fSubDetector2() const
Definition: JetIDHelper.h:50
double fLSbad() const
Definition: JetIDHelper.h:60
double fSubDetector4() const
Definition: JetIDHelper.h:52
reco::helper::JetIDHelper helper_
Definition: JetIDProducer.h:54
void calculate(const edm::Event &event, const edm::EventSetup &setup, const reco::CaloJet &jet, const int iDbg=0)
Definition: JetIDHelper.cc:93
void calculate(const edm::Event &event, const edm::EventSetup &isetup, const reco::Jet &jet, const int iDbg=0)
int iEvent
Definition: GenABIO.cc:224
double fHFOOT() const
Definition: JetIDHelper.h:61
double approximatefHPD() const
Definition: JetIDHelper.h:67
reco::helper::JetMuonHitsIDHelper muHelper_
Definition: JetIDProducer.h:55
double fShort() const
Definition: JetIDHelper.h:59
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
edm::EDGetTokenT< edm::View< reco::CaloJet > > input_jet_token_
Definition: JetIDProducer.h:57
double fSubDetector3() const
Definition: JetIDHelper.h:51
double restrictedEMF() const
Definition: JetIDHelper.h:63
double fLong() const
Definition: JetIDHelper.h:58
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

◆ helper_

reco::helper::JetIDHelper JetIDProducer::helper_
private

Definition at line 54 of file JetIDProducer.h.

Referenced by produce().

◆ input_jet_token_

edm::EDGetTokenT<edm::View<reco::CaloJet> > JetIDProducer::input_jet_token_
private

Definition at line 57 of file JetIDProducer.h.

Referenced by JetIDProducer(), and produce().

◆ muHelper_

reco::helper::JetMuonHitsIDHelper JetIDProducer::muHelper_
private

Definition at line 55 of file JetIDProducer.h.

Referenced by produce().

◆ src_

edm::InputTag JetIDProducer::src_
private

Definition at line 53 of file JetIDProducer.h.

Referenced by JetIDProducer().