Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:18:32

0001 #include "RecoJets/JetProducers/plugins/JetIDProducer.h"
0002 #include "DataFormats/JetReco/interface/JetID.h"
0003 
0004 #include <vector>
0005 
0006 //
0007 // constants, enums and typedefs
0008 //
0009 
0010 //
0011 // static data member definitions
0012 //
0013 
0014 //
0015 // constructors and destructor
0016 //
0017 JetIDProducer::JetIDProducer(const edm::ParameterSet& iConfig)
0018     : src_(iConfig.getParameter<edm::InputTag>("src")),
0019       helper_(iConfig, consumesCollector()),
0020       muHelper_(iConfig, consumesCollector()) {
0021   produces<reco::JetIDValueMap>();
0022 
0023   input_jet_token_ = consumes<edm::View<reco::CaloJet> >(src_);
0024 }
0025 
0026 JetIDProducer::~JetIDProducer() {}
0027 
0028 //
0029 // member functions
0030 //
0031 
0032 // ------------ method called to produce the data  ------------
0033 void JetIDProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0034   // get the input jets
0035   edm::Handle<edm::View<reco::CaloJet> > h_jets;
0036   iEvent.getByToken(input_jet_token_, h_jets);
0037 
0038   // allocate the jet--->jetid value map
0039   auto jetIdValueMap = std::make_unique<reco::JetIDValueMap>();
0040   // instantiate the filler with the map
0041   reco::JetIDValueMap::Filler filler(*jetIdValueMap);
0042 
0043   // allocate the vector of ids
0044   size_t njets = h_jets->size();
0045   std::vector<reco::JetID> ids(njets);
0046 
0047   // loop over the jets
0048   for (edm::View<reco::CaloJet>::const_iterator jetsBegin = h_jets->begin(), jetsEnd = h_jets->end(), ijet = jetsBegin;
0049        ijet != jetsEnd;
0050        ++ijet) {
0051     // get the id from each jet
0052     helper_.calculate(iEvent, iSetup, *ijet);
0053 
0054     muHelper_.calculate(iEvent, iSetup, *ijet);
0055 
0056     ids[ijet - jetsBegin].fHPD = helper_.fHPD();
0057     ids[ijet - jetsBegin].fRBX = helper_.fRBX();
0058     ids[ijet - jetsBegin].n90Hits = helper_.n90Hits();
0059     ids[ijet - jetsBegin].fSubDetector1 = helper_.fSubDetector1();
0060     ids[ijet - jetsBegin].fSubDetector2 = helper_.fSubDetector2();
0061     ids[ijet - jetsBegin].fSubDetector3 = helper_.fSubDetector3();
0062     ids[ijet - jetsBegin].fSubDetector4 = helper_.fSubDetector4();
0063     ids[ijet - jetsBegin].restrictedEMF = helper_.restrictedEMF();
0064     ids[ijet - jetsBegin].nHCALTowers = helper_.nHCALTowers();
0065     ids[ijet - jetsBegin].nECALTowers = helper_.nECALTowers();
0066     ids[ijet - jetsBegin].approximatefHPD = helper_.approximatefHPD();
0067     ids[ijet - jetsBegin].approximatefRBX = helper_.approximatefRBX();
0068     ids[ijet - jetsBegin].hitsInN90 = helper_.hitsInN90();
0069 
0070     ids[ijet - jetsBegin].numberOfHits2RPC = muHelper_.numberOfHits2RPC();
0071     ids[ijet - jetsBegin].numberOfHits3RPC = muHelper_.numberOfHits3RPC();
0072     ids[ijet - jetsBegin].numberOfHitsRPC = muHelper_.numberOfHitsRPC();
0073 
0074     ids[ijet - jetsBegin].fEB = helper_.fEB();
0075     ids[ijet - jetsBegin].fEE = helper_.fEE();
0076     ids[ijet - jetsBegin].fHB = helper_.fHB();
0077     ids[ijet - jetsBegin].fHE = helper_.fHE();
0078     ids[ijet - jetsBegin].fHO = helper_.fHO();
0079     ids[ijet - jetsBegin].fLong = helper_.fLong();
0080     ids[ijet - jetsBegin].fShort = helper_.fShort();
0081     ids[ijet - jetsBegin].fLS = helper_.fLSbad();
0082     ids[ijet - jetsBegin].fHFOOT = helper_.fHFOOT();
0083   }
0084 
0085   // set up the map
0086   filler.insert(h_jets, ids.begin(), ids.end());
0087 
0088   // fill the vals
0089   filler.fill();
0090 
0091   // write map to the event
0092   iEvent.put(std::move(jetIdValueMap));
0093 }
0094 
0095 //define this as a plug-in
0096 DEFINE_FWK_MODULE(JetIDProducer);