Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:23:47

0001 
0002 #include "FWCore/Framework/interface/stream/EDProducer.h"
0003 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0004 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0005 #include "DataFormats/Common/interface/RefVector.h"
0006 
0007 #include "CommonTools/UtilAlgos/interface/StringCutObjectSelector.h"
0008 #include "CommonTools/UtilAlgos/interface/SingleObjectSelector.h"
0009 #include "CommonTools/UtilAlgos/interface/ObjectSelector.h"
0010 #include "CommonTools/UtilAlgos/interface/SingleElementCollectionSelector.h"
0011 
0012 #include "DataFormats/JetReco/interface/JPTJet.h"
0013 #include "DataFormats/JetReco/interface/JPTJetCollection.h"
0014 #include "DataFormats/JetReco/interface/CaloJet.h"
0015 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
0016 #include "DataFormats/Candidate/interface/Candidate.h"
0017 #include "DataFormats/Math/interface/deltaR.h"
0018 #include <vector>
0019 
0020 class JPTJetSlimmer : public edm::stream::EDProducer<> {
0021 public:
0022   JPTJetSlimmer(edm::ParameterSet const& params)
0023       : srcToken_(consumes<edm::View<reco::JPTJet> >(params.getParameter<edm::InputTag>("src"))),
0024         srcCaloToken_(consumes<edm::View<reco::CaloJet> >(params.getParameter<edm::InputTag>("srcCalo"))),
0025         cut_(params.getParameter<std::string>("cut")),
0026         selector_(cut_) {
0027     produces<reco::JPTJetCollection>();
0028     produces<reco::CaloJetCollection>();
0029   }
0030 
0031   ~JPTJetSlimmer() override {}
0032 
0033   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0034 
0035   void produce(edm::Event& iEvent, const edm::EventSetup& iSetup) override {
0036     auto jptJets = std::make_unique<reco::JPTJetCollection>();
0037     auto caloJets = std::make_unique<reco::CaloJetCollection>();
0038 
0039     edm::RefProd<reco::CaloJetCollection> pOut1RefProd = iEvent.getRefBeforePut<reco::CaloJetCollection>();
0040     edm::Ref<reco::CaloJetCollection>::key_type idxCaloJet = 0;
0041 
0042     auto const& h_calojets = iEvent.get(srcCaloToken_);
0043     auto const& h_jets = iEvent.get(srcToken_);
0044 
0045     for (auto const& ijet : h_jets) {
0046       if (selector_(ijet)) {
0047         // Add specific : only reference to CaloJet collection. It is necessary for
0048         // recalibration JPTJet at MiniAod.
0049         const edm::RefToBase<reco::Jet>& rawcalojet = ijet.getCaloJetRef();
0050         int icalo = -1;
0051         int i = 0;
0052         for (auto const& icjet : h_calojets) {
0053           double dr2 = deltaR2(icjet, *rawcalojet);
0054           if (dr2 <= 0.001) {
0055             icalo = i;
0056           }
0057           i++;
0058         }
0059         reco::JPTJet::Specific tmp_specific;
0060         if (icalo < 0) {
0061           // Add reference to the created CaloJet collection
0062           reco::CaloJetRef myjet(pOut1RefProd, idxCaloJet++);
0063           tmp_specific.theCaloJetRef = edm::RefToBase<reco::Jet>(myjet);
0064           reco::CaloJet const* rawcalojetc = dynamic_cast<reco::CaloJet const*>(&*rawcalojet);
0065           caloJets->push_back(*rawcalojetc);
0066           const reco::Candidate::Point& orivtx = ijet.vertex();
0067           reco::JPTJet newJPTJet(ijet.p4(), orivtx, tmp_specific, ijet.getJetConstituents());
0068           float jetArea = ijet.jetArea();
0069           newJPTJet.setJetArea(fabs(jetArea));
0070           jptJets->push_back(newJPTJet);
0071         }
0072       }
0073     }
0074     iEvent.put(std::move(caloJets));
0075     iEvent.put(std::move(jptJets));
0076   }
0077 
0078 protected:
0079   const edm::EDGetTokenT<edm::View<reco::JPTJet> > srcToken_;
0080   const edm::EDGetTokenT<edm::View<reco::CaloJet> > srcCaloToken_;
0081   const std::string cut_;
0082   const StringCutObjectSelector<reco::Jet> selector_;
0083 };
0084 void JPTJetSlimmer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0085   // slimmedJPTJets
0086   edm::ParameterSetDescription desc;
0087   desc.add<edm::InputTag>("src", edm::InputTag("JetPlusTrackZSPCorJetAntiKt4"));
0088   desc.add<edm::InputTag>("srcCalo", edm::InputTag("slimmedCaloJets"));
0089   desc.add<std::string>("cut", "pt>20");
0090   descriptions.add("slimmedJPTJets", desc);
0091 }
0092 
0093 #include "FWCore/Framework/interface/MakerMacros.h"
0094 DEFINE_FWK_MODULE(JPTJetSlimmer);