Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:01:05

0001 // Weight for neutral particles based on distance with charged
0002 //
0003 // Original Author:  Michail Bachtis,40 1-B08,+41227678176,
0004 //         Created:  Mon Dec  9 13:18:05 CET 2013
0005 //
0006 // edited by Pavel Jez
0007 //
0008 
0009 #include "DataFormats/Math/interface/deltaR.h"
0010 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0011 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
0012 #include "FWCore/Framework/interface/global/EDProducer.h"
0013 #include "FWCore/Framework/interface/Event.h"
0014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0016 
0017 #include <memory>
0018 
0019 class DeltaBetaWeights : public edm::global::EDProducer<> {
0020 public:
0021   explicit DeltaBetaWeights(const edm::ParameterSet&);
0022 
0023 private:
0024   void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0025   // ----------member data ---------------------------
0026   edm::InputTag src_;
0027   edm::InputTag pfCharged_;
0028   edm::InputTag pfPU_;
0029 
0030   edm::EDGetTokenT<edm::View<reco::Candidate> > pfCharged_token;
0031   edm::EDGetTokenT<edm::View<reco::Candidate> > pfPU_token;
0032   edm::EDGetTokenT<edm::View<reco::Candidate> > src_token;
0033 };
0034 
0035 DeltaBetaWeights::DeltaBetaWeights(const edm::ParameterSet& iConfig)
0036     : src_(iConfig.getParameter<edm::InputTag>("src")),
0037       pfCharged_(iConfig.getParameter<edm::InputTag>("chargedFromPV")),
0038       pfPU_(iConfig.getParameter<edm::InputTag>("chargedFromPU")) {
0039   produces<reco::PFCandidateCollection>();
0040 
0041   pfCharged_token = consumes<edm::View<reco::Candidate> >(pfCharged_);
0042   pfPU_token = consumes<edm::View<reco::Candidate> >(pfPU_);
0043   src_token = consumes<edm::View<reco::Candidate> >(src_);
0044 
0045   // pfCharged_token = consumes<reco::PFCandidateCollection>(pfCharged_);
0046   // pfPU_token = consumes<reco::PFCandidateCollection>(pfPU_);
0047   // src_token = consumes<reco::PFCandidateCollection>(src_);
0048 }
0049 
0050 #include "FWCore/Framework/interface/MakerMacros.h"
0051 DEFINE_FWK_MODULE(DeltaBetaWeights);
0052 
0053 // ------------ method called to produce the data  ------------
0054 void DeltaBetaWeights::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0055   using namespace edm;
0056 
0057   edm::View<reco::Candidate> const pfCharged = iEvent.get(pfCharged_token);
0058   edm::View<reco::Candidate> const& pfPU = iEvent.get(pfPU_token);
0059   edm::View<reco::Candidate> const& src = iEvent.get(src_token);
0060 
0061   double sumNPU = .0;
0062   double sumPU = .0;
0063 
0064   std::unique_ptr<reco::PFCandidateCollection> out(new reco::PFCandidateCollection);
0065 
0066   out->reserve(src.size());
0067   for (const reco::Candidate& cand : src) {
0068     if (cand.charge() != 0) {
0069       // this part of code should be executed only if input collection is not entirely composed of neutral candidates, i.e. never by default
0070       edm::LogWarning("DeltaBetaWeights")
0071           << "Trying to reweight charged particle... saving it to output collection without any change";
0072       out->emplace_back(cand.charge(), cand.p4(), reco::PFCandidate::ParticleType::X);
0073       (out->back()).setParticleType((out->back()).translatePdgIdToType(cand.pdgId()));
0074       continue;
0075     }
0076 
0077     sumNPU = 1.0;
0078     sumPU = 1.0;
0079     double eta = cand.eta();
0080     double phi = cand.phi();
0081     for (const reco::Candidate& chCand : pfCharged) {
0082       double sum = (chCand.pt() * chCand.pt()) / (deltaR2(eta, phi, chCand.eta(), chCand.phi()));
0083       if (sum > 1.0)
0084         sumNPU *= sum;
0085     }
0086     sumNPU = 0.5 * log(sumNPU);
0087 
0088     for (const reco::Candidate& puCand : pfPU) {
0089       double sum = (puCand.pt() * puCand.pt()) / (deltaR2(eta, phi, puCand.eta(), puCand.phi()));
0090       if (sum > 1.0)
0091         sumPU *= sum;
0092     }
0093     sumPU = 0.5 * log(sumPU);
0094 
0095     reco::PFCandidate neutral = reco::PFCandidate(cand.charge(), cand.p4(), reco::PFCandidate::ParticleType::X);
0096     neutral.setParticleType(neutral.translatePdgIdToType(cand.pdgId()));
0097     if (sumNPU + sumPU > 0)
0098       neutral.setP4(((sumNPU) / (sumNPU + sumPU)) * neutral.p4());
0099     out->push_back(neutral);
0100   }
0101 
0102   iEvent.put(std::move(out));
0103 }