CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoJets/JetProducers/plugins/JetIDProducer.cc

Go to the documentation of this file.
00001 #include "RecoJets/JetProducers/plugins/JetIDProducer.h"
00002 #include "DataFormats/JetReco/interface/JetID.h"
00003 
00004 #include <vector>
00005 
00006 //
00007 // constants, enums and typedefs
00008 //
00009 
00010 
00011 //
00012 // static data member definitions
00013 //
00014 
00015 //
00016 // constructors and destructor
00017 //
00018 JetIDProducer::JetIDProducer(const edm::ParameterSet& iConfig) :
00019   src_       ( iConfig.getParameter<edm::InputTag>("src") ),
00020   helper_    ( iConfig ),
00021   muHelper_  ( iConfig )
00022 {
00023   produces< reco::JetIDValueMap >();
00024 }
00025 
00026 
00027 JetIDProducer::~JetIDProducer()
00028 {
00029 }
00030 
00031 
00032 //
00033 // member functions
00034 //
00035 
00036 // ------------ method called to produce the data  ------------
00037 void
00038 JetIDProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00039 {
00040 
00041   // get the input jets
00042   edm::Handle< edm::View<reco::CaloJet> > h_jets;
00043   iEvent.getByLabel( src_, h_jets );
00044 
00045   // allocate the jet--->jetid value map
00046   std::auto_ptr<reco::JetIDValueMap> jetIdValueMap( new reco::JetIDValueMap );
00047   // instantiate the filler with the map
00048   reco::JetIDValueMap::Filler filler(*jetIdValueMap);
00049   
00050   // allocate the vector of ids
00051   size_t njets = h_jets->size();
00052   std::vector<reco::JetID>  ids (njets);
00053    
00054   // loop over the jets
00055   for ( edm::View<reco::CaloJet>::const_iterator jetsBegin = h_jets->begin(),
00056           jetsEnd = h_jets->end(),
00057           ijet = jetsBegin;
00058         ijet != jetsEnd; ++ijet ) {
00059 
00060     // get the id from each jet
00061     helper_.calculate( iEvent, *ijet );
00062 
00063     muHelper_.calculate( iEvent, iSetup, *ijet );
00064 
00065     ids[ijet-jetsBegin].fHPD               =  helper_.fHPD();
00066     ids[ijet-jetsBegin].fRBX               =  helper_.fRBX();
00067     ids[ijet-jetsBegin].n90Hits            =  helper_.n90Hits();
00068     ids[ijet-jetsBegin].fSubDetector1      =  helper_.fSubDetector1();
00069     ids[ijet-jetsBegin].fSubDetector2      =  helper_.fSubDetector2();
00070     ids[ijet-jetsBegin].fSubDetector3      =  helper_.fSubDetector3();
00071     ids[ijet-jetsBegin].fSubDetector4      =  helper_.fSubDetector4();
00072     ids[ijet-jetsBegin].restrictedEMF      =  helper_.restrictedEMF();
00073     ids[ijet-jetsBegin].nHCALTowers        =  helper_.nHCALTowers();
00074     ids[ijet-jetsBegin].nECALTowers        =  helper_.nECALTowers();
00075     ids[ijet-jetsBegin].approximatefHPD    =  helper_.approximatefHPD();
00076     ids[ijet-jetsBegin].approximatefRBX    =  helper_.approximatefRBX();
00077     ids[ijet-jetsBegin].hitsInN90          =  helper_.hitsInN90();    
00078 
00079     ids[ijet-jetsBegin].numberOfHits2RPC   = muHelper_.numberOfHits2RPC();
00080     ids[ijet-jetsBegin].numberOfHits3RPC   = muHelper_.numberOfHits3RPC();
00081     ids[ijet-jetsBegin].numberOfHitsRPC    = muHelper_.numberOfHitsRPC();
00082     
00083     ids[ijet-jetsBegin].fEB     = helper_.fEB   ();
00084     ids[ijet-jetsBegin].fEE     = helper_.fEE   ();
00085     ids[ijet-jetsBegin].fHB     = helper_.fHB   (); 
00086     ids[ijet-jetsBegin].fHE     = helper_.fHE   (); 
00087     ids[ijet-jetsBegin].fHO     = helper_.fHO   (); 
00088     ids[ijet-jetsBegin].fLong   = helper_.fLong ();
00089     ids[ijet-jetsBegin].fShort  = helper_.fShort();
00090     ids[ijet-jetsBegin].fLS     = helper_.fLSbad   ();
00091     ids[ijet-jetsBegin].fHFOOT  = helper_.fHFOOT();
00092   }
00093   
00094   // set up the map
00095   filler.insert( h_jets, ids.begin(), ids.end() );
00096  
00097   // fill the vals
00098   filler.fill();
00099 
00100   // write map to the event
00101   iEvent.put( jetIdValueMap );
00102 }
00103 
00104 // ------------ method called once each job just before starting event loop  ------------
00105 void 
00106 JetIDProducer::beginJob()
00107 {
00108 }
00109 
00110 // ------------ method called once each job just after ending the event loop  ------------
00111 void 
00112 JetIDProducer::endJob() {
00113 }
00114 
00115 //define this as a plug-in
00116 DEFINE_FWK_MODULE(JetIDProducer);