CMS 3D CMS Logo

DeltaBetaWeights.cc
Go to the documentation of this file.
1 // Weight for neutral particles based on distance with charged
2 //
3 // Original Author: Michail Bachtis,40 1-B08,+41227678176,
4 // Created: Mon Dec 9 13:18:05 CET 2013
5 //
6 // edited by Pavel Jez
7 //
8 
16 
17 #include <memory>
18 
20 public:
21  explicit DeltaBetaWeights(const edm::ParameterSet&);
22 
23 private:
24  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
25  // ----------member data ---------------------------
29 
33 };
34 
36  : src_(iConfig.getParameter<edm::InputTag>("src")),
37  pfCharged_(iConfig.getParameter<edm::InputTag>("chargedFromPV")),
38  pfPU_(iConfig.getParameter<edm::InputTag>("chargedFromPU")) {
39  produces<reco::PFCandidateCollection>();
40 
41  pfCharged_token = consumes<edm::View<reco::Candidate> >(pfCharged_);
42  pfPU_token = consumes<edm::View<reco::Candidate> >(pfPU_);
43  src_token = consumes<edm::View<reco::Candidate> >(src_);
44 
45  // pfCharged_token = consumes<reco::PFCandidateCollection>(pfCharged_);
46  // pfPU_token = consumes<reco::PFCandidateCollection>(pfPU_);
47  // src_token = consumes<reco::PFCandidateCollection>(src_);
48 }
49 
52 
53 // ------------ method called to produce the data ------------
55  using namespace edm;
56 
57  edm::View<reco::Candidate> const pfCharged = iEvent.get(pfCharged_token);
58  edm::View<reco::Candidate> const& pfPU = iEvent.get(pfPU_token);
60 
61  double sumNPU = .0;
62  double sumPU = .0;
63 
64  std::unique_ptr<reco::PFCandidateCollection> out(new reco::PFCandidateCollection);
65 
66  out->reserve(src.size());
67  for (const reco::Candidate& cand : src) {
68  if (cand.charge() != 0) {
69  // this part of code should be executed only if input collection is not entirely composed of neutral candidates, i.e. never by default
70  edm::LogWarning("DeltaBetaWeights")
71  << "Trying to reweight charged particle... saving it to output collection without any change";
72  out->emplace_back(cand.charge(), cand.p4(), reco::PFCandidate::ParticleType::X);
73  (out->back()).setParticleType((out->back()).translatePdgIdToType(cand.pdgId()));
74  continue;
75  }
76 
77  sumNPU = 1.0;
78  sumPU = 1.0;
79  double eta = cand.eta();
80  double phi = cand.phi();
81  for (const reco::Candidate& chCand : pfCharged) {
82  double sum = (chCand.pt() * chCand.pt()) / (deltaR2(eta, phi, chCand.eta(), chCand.phi()));
83  if (sum > 1.0)
84  sumNPU *= sum;
85  }
86  sumNPU = 0.5 * log(sumNPU);
87 
88  for (const reco::Candidate& puCand : pfPU) {
89  double sum = (puCand.pt() * puCand.pt()) / (deltaR2(eta, phi, puCand.eta(), puCand.phi()));
90  if (sum > 1.0)
91  sumPU *= sum;
92  }
93  sumPU = 0.5 * log(sumPU);
94 
96  neutral.setParticleType(neutral.translatePdgIdToType(cand.pdgId()));
97  if (sumNPU + sumPU > 0)
98  neutral.setP4(((sumNPU) / (sumNPU + sumPU)) * neutral.p4());
99  out->push_back(neutral);
100  }
101 
102  iEvent.put(std::move(out));
103 }
edm::EDGetTokenT< edm::View< reco::Candidate > > pfPU_token
edm::EDGetTokenT< edm::View< reco::Candidate > > pfCharged_token
#define X(str)
Definition: MuonsGrabber.cc:38
const LorentzVector & p4() const final
four-momentum Lorentz vector
int iEvent
Definition: GenABIO.cc:224
void setParticleType(ParticleType type)
set Particle Type
Definition: PFCandidate.cc:278
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
edm::InputTag pfCharged_
edm::InputTag pfPU_
edm::EDGetTokenT< edm::View< reco::Candidate > > src_token
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
ParticleType translatePdgIdToType(int pdgid) const
Definition: PFCandidate.cc:231
HLT enums.
Log< level::Warning, false > LogWarning
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
edm::InputTag src_
void setP4(const LorentzVector &p4) final
set 4-momentum
def move(src, dest)
Definition: eostools.py:511
DeltaBetaWeights(const edm::ParameterSet &)