CMS 3D CMS Logo

JetIDProducer.cc
Go to the documentation of this file.
3 
4 #include <vector>
5 
6 //
7 // constants, enums and typedefs
8 //
9 
10 
11 //
12 // static data member definitions
13 //
14 
15 //
16 // constructors and destructor
17 //
19  src_ ( iConfig.getParameter<edm::InputTag>("src") ),
20  helper_ ( iConfig, consumesCollector() ),
21  muHelper_ ( iConfig, consumesCollector() )
22 {
23  produces< reco::JetIDValueMap >();
24 
25  input_jet_token_ = consumes<edm::View<reco::CaloJet> >(src_);
26 
27 }
28 
29 
31 {
32 }
33 
34 
35 //
36 // member functions
37 //
38 
39 // ------------ method called to produce the data ------------
40 void
42 {
43 
44  // get the input jets
46  iEvent.getByToken( input_jet_token_, h_jets );
47 
48  // allocate the jet--->jetid value map
49  auto jetIdValueMap = std::make_unique<reco::JetIDValueMap>();
50  // instantiate the filler with the map
51  reco::JetIDValueMap::Filler filler(*jetIdValueMap);
52 
53  // allocate the vector of ids
54  size_t njets = h_jets->size();
55  std::vector<reco::JetID> ids (njets);
56 
57  // loop over the jets
58  for ( edm::View<reco::CaloJet>::const_iterator jetsBegin = h_jets->begin(),
59  jetsEnd = h_jets->end(),
60  ijet = jetsBegin;
61  ijet != jetsEnd; ++ijet ) {
62 
63  // get the id from each jet
64  helper_.calculate( iEvent, iSetup, *ijet );
65 
66  muHelper_.calculate( iEvent, iSetup, *ijet );
67 
68  ids[ijet-jetsBegin].fHPD = helper_.fHPD();
69  ids[ijet-jetsBegin].fRBX = helper_.fRBX();
70  ids[ijet-jetsBegin].n90Hits = helper_.n90Hits();
71  ids[ijet-jetsBegin].fSubDetector1 = helper_.fSubDetector1();
72  ids[ijet-jetsBegin].fSubDetector2 = helper_.fSubDetector2();
73  ids[ijet-jetsBegin].fSubDetector3 = helper_.fSubDetector3();
74  ids[ijet-jetsBegin].fSubDetector4 = helper_.fSubDetector4();
75  ids[ijet-jetsBegin].restrictedEMF = helper_.restrictedEMF();
76  ids[ijet-jetsBegin].nHCALTowers = helper_.nHCALTowers();
77  ids[ijet-jetsBegin].nECALTowers = helper_.nECALTowers();
78  ids[ijet-jetsBegin].approximatefHPD = helper_.approximatefHPD();
79  ids[ijet-jetsBegin].approximatefRBX = helper_.approximatefRBX();
80  ids[ijet-jetsBegin].hitsInN90 = helper_.hitsInN90();
81 
82  ids[ijet-jetsBegin].numberOfHits2RPC = muHelper_.numberOfHits2RPC();
83  ids[ijet-jetsBegin].numberOfHits3RPC = muHelper_.numberOfHits3RPC();
84  ids[ijet-jetsBegin].numberOfHitsRPC = muHelper_.numberOfHitsRPC();
85 
86  ids[ijet-jetsBegin].fEB = helper_.fEB ();
87  ids[ijet-jetsBegin].fEE = helper_.fEE ();
88  ids[ijet-jetsBegin].fHB = helper_.fHB ();
89  ids[ijet-jetsBegin].fHE = helper_.fHE ();
90  ids[ijet-jetsBegin].fHO = helper_.fHO ();
91  ids[ijet-jetsBegin].fLong = helper_.fLong ();
92  ids[ijet-jetsBegin].fShort = helper_.fShort();
93  ids[ijet-jetsBegin].fLS = helper_.fLSbad ();
94  ids[ijet-jetsBegin].fHFOOT = helper_.fHFOOT();
95  }
96 
97  // set up the map
98  filler.insert( h_jets, ids.begin(), ids.end() );
99 
100  // fill the vals
101  filler.fill();
102 
103  // write map to the event
104  iEvent.put(std::move(jetIdValueMap));
105 }
106 
107 //define this as a plug-in
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:136
double fHE() const
Definition: JetIDHelper.h:51
double fHB() const
Definition: JetIDHelper.h:50
double fLSbad() const
Definition: JetIDHelper.h:55
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
int nHCALTowers() const
Definition: JetIDHelper.h:59
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
double approximatefRBX() const
Definition: JetIDHelper.h:63
double fLong() const
Definition: JetIDHelper.h:53
int nECALTowers() const
Definition: JetIDHelper.h:60
double fHPD() const
Definition: JetIDHelper.h:41
double fSubDetector4() const
Definition: JetIDHelper.h:47
reco::helper::JetIDHelper helper_
Definition: JetIDProducer.h:56
double restrictedEMF() const
Definition: JetIDHelper.h:58
void calculate(const edm::Event &event, const edm::EventSetup &setup, const reco::CaloJet &jet, const int iDbg=0)
Definition: JetIDHelper.cc:100
double fSubDetector1() const
Definition: JetIDHelper.h:44
~JetIDProducer() override
void calculate(const edm::Event &event, const edm::EventSetup &isetup, const reco::Jet &jet, const int iDbg=0)
edm::InputTag src_
Definition: JetIDProducer.h:55
int iEvent
Definition: GenABIO.cc:230
double approximatefHPD() const
Definition: JetIDHelper.h:62
double fEB() const
Definition: JetIDHelper.h:48
JetIDProducer(const edm::ParameterSet &)
double fEE() const
Definition: JetIDHelper.h:49
double fSubDetector3() const
Definition: JetIDHelper.h:46
double fShort() const
Definition: JetIDHelper.h:54
reco::helper::JetMuonHitsIDHelper muHelper_
Definition: JetIDProducer.h:57
HLT enums.
double fSubDetector2() const
Definition: JetIDHelper.h:45
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
void produce(edm::Event &, const edm::EventSetup &) override
edm::EDGetTokenT< edm::View< reco::CaloJet > > input_jet_token_
Definition: JetIDProducer.h:59
double fHO() const
Definition: JetIDHelper.h:52
double fHFOOT() const
Definition: JetIDHelper.h:56
def move(src, dest)
Definition: eostools.py:510
double fRBX() const
Definition: JetIDHelper.h:42