Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:31

0001 // -*- C++ -*-
0002 //
0003 // Package:    BasicToPFJet
0004 // Class:      BasicToPFJet
0005 //
0006 /**\class BasicToPFJet BasicToPFJet.cc UserCode/BasicToPFJet/plugins/BasicToPFJet.cc
0007 
0008  Description: converts reco::BasicJets to reco::PFJets and adds the new PFJetCollection to the event. Originally designed
0009               to be a work around for a way to store reco::BasicJets at HLT level
0010 
0011  Implementation:
0012      [Notes on implementation]
0013 */
0014 //
0015 // Original Author:  clint richardson
0016 //         Created:  Thu, 6 Mar 2014 12:00:00 GMT
0017 //
0018 //
0019 // system include files
0020 #include <memory>
0021 #include <vector>
0022 #include <sstream>
0023 
0024 // user include files
0025 #include "FWCore/PluginManager/interface/ModuleDef.h"
0026 #include "FWCore/Framework/interface/MakerMacros.h"
0027 #include "FWCore/Framework/interface/Event.h"
0028 #include "DataFormats/BTauReco/interface/CATopJetTagInfo.h"
0029 #include "DataFormats/JetReco/interface/PFJetCollection.h"
0030 
0031 //include header file
0032 #include "RecoJets/JetProducers/plugins/BasicToPFJet.h"
0033 
0034 BasicToPFJet::BasicToPFJet(const edm::ParameterSet& PSet)
0035     : src_(PSet.getParameter<edm::InputTag>("src")), inputToken_(consumes<reco::BasicJetCollection>(src_)) {
0036   produces<reco::PFJetCollection>();
0037 }
0038 
0039 void BasicToPFJet::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0040   edm::ParameterSetDescription desc;
0041   desc.add<edm::InputTag>("src", edm::InputTag(""));
0042   descriptions.add("BasicToPFJet", desc);
0043 }
0044 
0045 void BasicToPFJet::produce(edm::StreamID, edm::Event& Event, const edm::EventSetup& EventSetup) const {
0046   //first get the basic jet collection
0047   edm::Handle<reco::BasicJetCollection> BasicJetColl;
0048   Event.getByToken(inputToken_, BasicJetColl);
0049 
0050   //now make the new pf jet collection
0051   auto PFJetColl = std::make_unique<reco::PFJetCollection>();
0052   //reco::PFJetCollection* PFJetColl = new reco::PFJetCollection;
0053   //make the 'specific'
0054   reco::PFJet::Specific specific;
0055 
0056   //now get iterator
0057   reco::BasicJetCollection::const_iterator i = BasicJetColl->begin();
0058 
0059   //loop over basic jets and convert them to pfjets
0060   for (; i != BasicJetColl->end(); i++) {
0061     reco::PFJet pfjet(i->p4(), i->vertex(), specific);
0062     PFJetColl->push_back(pfjet);
0063   }
0064 
0065   //std::unique_ptr<reco::PFJetCollection> selectedPFJets(PFJetColl);
0066   Event.put(std::move(PFJetColl));
0067 }
0068 
0069 //define as plug-in for the framework
0070 DEFINE_FWK_MODULE(BasicToPFJet);