1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
|
// Weight for neutral particles based on distance with charged
//
// Original Author: Michail Bachtis,40 1-B08,+41227678176,
// Created: Mon Dec 9 13:18:05 CET 2013
//
// edited by Pavel Jez
//
#include "DataFormats/Math/interface/deltaR.h"
#include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
#include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
#include "FWCore/Framework/interface/global/EDProducer.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include <memory>
class DeltaBetaWeights : public edm::global::EDProducer<> {
public:
explicit DeltaBetaWeights(const edm::ParameterSet&);
private:
void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
// ----------member data ---------------------------
edm::InputTag src_;
edm::InputTag pfCharged_;
edm::InputTag pfPU_;
edm::EDGetTokenT<edm::View<reco::Candidate> > pfCharged_token;
edm::EDGetTokenT<edm::View<reco::Candidate> > pfPU_token;
edm::EDGetTokenT<edm::View<reco::Candidate> > src_token;
};
DeltaBetaWeights::DeltaBetaWeights(const edm::ParameterSet& iConfig)
: src_(iConfig.getParameter<edm::InputTag>("src")),
pfCharged_(iConfig.getParameter<edm::InputTag>("chargedFromPV")),
pfPU_(iConfig.getParameter<edm::InputTag>("chargedFromPU")) {
produces<reco::PFCandidateCollection>();
pfCharged_token = consumes<edm::View<reco::Candidate> >(pfCharged_);
pfPU_token = consumes<edm::View<reco::Candidate> >(pfPU_);
src_token = consumes<edm::View<reco::Candidate> >(src_);
// pfCharged_token = consumes<reco::PFCandidateCollection>(pfCharged_);
// pfPU_token = consumes<reco::PFCandidateCollection>(pfPU_);
// src_token = consumes<reco::PFCandidateCollection>(src_);
}
#include "FWCore/Framework/interface/MakerMacros.h"
DEFINE_FWK_MODULE(DeltaBetaWeights);
// ------------ method called to produce the data ------------
void DeltaBetaWeights::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
using namespace edm;
edm::View<reco::Candidate> const pfCharged = iEvent.get(pfCharged_token);
edm::View<reco::Candidate> const& pfPU = iEvent.get(pfPU_token);
edm::View<reco::Candidate> const& src = iEvent.get(src_token);
double sumNPU = .0;
double sumPU = .0;
std::unique_ptr<reco::PFCandidateCollection> out(new reco::PFCandidateCollection);
out->reserve(src.size());
for (const reco::Candidate& cand : src) {
if (cand.charge() != 0) {
// this part of code should be executed only if input collection is not entirely composed of neutral candidates, i.e. never by default
edm::LogWarning("DeltaBetaWeights")
<< "Trying to reweight charged particle... saving it to output collection without any change";
out->emplace_back(cand.charge(), cand.p4(), reco::PFCandidate::ParticleType::X);
(out->back()).setParticleType((out->back()).translatePdgIdToType(cand.pdgId()));
continue;
}
sumNPU = 1.0;
sumPU = 1.0;
double eta = cand.eta();
double phi = cand.phi();
for (const reco::Candidate& chCand : pfCharged) {
double sum = (chCand.pt() * chCand.pt()) / (deltaR2(eta, phi, chCand.eta(), chCand.phi()));
if (sum > 1.0)
sumNPU *= sum;
}
sumNPU = 0.5 * log(sumNPU);
for (const reco::Candidate& puCand : pfPU) {
double sum = (puCand.pt() * puCand.pt()) / (deltaR2(eta, phi, puCand.eta(), puCand.phi()));
if (sum > 1.0)
sumPU *= sum;
}
sumPU = 0.5 * log(sumPU);
reco::PFCandidate neutral = reco::PFCandidate(cand.charge(), cand.p4(), reco::PFCandidate::ParticleType::X);
neutral.setParticleType(neutral.translatePdgIdToType(cand.pdgId()));
if (sumNPU + sumPU > 0)
neutral.setP4(((sumNPU) / (sumNPU + sumPU)) * neutral.p4());
out->push_back(neutral);
}
iEvent.put(std::move(out));
}
|