CMS 3D CMS Logo

DeltaBetaWeights.cc
Go to the documentation of this file.
3 
5  : src_(iConfig.getParameter<edm::InputTag>("src")),
6  pfCharged_(iConfig.getParameter<edm::InputTag>("chargedFromPV")),
7  pfPU_(iConfig.getParameter<edm::InputTag>("chargedFromPU")) {
8  produces<reco::PFCandidateCollection>();
9 
10  pfCharged_token = consumes<edm::View<reco::Candidate> >(pfCharged_);
11  pfPU_token = consumes<edm::View<reco::Candidate> >(pfPU_);
12  src_token = consumes<edm::View<reco::Candidate> >(src_);
13 
14  // pfCharged_token = consumes<reco::PFCandidateCollection>(pfCharged_);
15  // pfPU_token = consumes<reco::PFCandidateCollection>(pfPU_);
16  // src_token = consumes<reco::PFCandidateCollection>(src_);
17 }
18 
20  // do anything here that needs to be done at desctruction time
21  // (e.g. close files, deallocate resources etc.)
22 }
23 
24 //
25 // member functions
26 //
27 
28 // ------------ method called to produce the data ------------
30  using namespace edm;
31 
35 
36  iEvent.getByToken(src_token, src);
37  iEvent.getByToken(pfCharged_token, pfCharged);
38  iEvent.getByToken(pfPU_token, pfPU);
39 
40  double sumNPU = .0;
41  double sumPU = .0;
42 
43  std::unique_ptr<reco::PFCandidateCollection> out(new reco::PFCandidateCollection);
44 
45  for (const reco::Candidate& cand : *src) {
46  if (cand.charge() != 0) {
47  // this part of code should be executed only if input collection is not entirely composed of neutral candidates, i.e. never by default
48  edm::LogWarning("DeltaBetaWeights")
49  << "Trying to reweight charged particle... saving it to output collection without any change";
50  out->emplace_back(cand.charge(), cand.p4(), reco::PFCandidate::ParticleType::X);
51  (out->back()).setParticleType((out->back()).translatePdgIdToType(cand.pdgId()));
52  continue;
53  }
54 
55  sumNPU = 1.0;
56  sumPU = 1.0;
57  double eta = cand.eta();
58  double phi = cand.phi();
59  for (const reco::Candidate& chCand : *pfCharged) {
60  double sum = (chCand.pt() * chCand.pt()) / (deltaR2(eta, phi, chCand.eta(), chCand.phi()));
61  if (sum > 1.0)
62  sumNPU *= sum;
63  }
64  sumNPU = 0.5 * log(sumNPU);
65 
66  for (const reco::Candidate& puCand : *pfPU) {
67  double sum = (puCand.pt() * puCand.pt()) / (deltaR2(eta, phi, puCand.eta(), puCand.phi()));
68  if (sum > 1.0)
69  sumPU *= sum;
70  }
71  sumPU = 0.5 * log(sumPU);
72 
74  neutral.setParticleType(neutral.translatePdgIdToType(cand.pdgId()));
75  if (sumNPU + sumPU > 0)
76  neutral.setP4(((sumNPU) / (sumNPU + sumPU)) * neutral.p4());
77  out->push_back(neutral);
78  }
79 
80  iEvent.put(std::move(out));
81 }
MessageLogger.h
X
#define X(str)
Definition: MuonsGrabber.cc:38
DeltaBetaWeights::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: DeltaBetaWeights.cc:29
edm
HLT enums.
Definition: AlignableModifier.h:19
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89285
DeltaBetaWeights::~DeltaBetaWeights
~DeltaBetaWeights() override
Definition: DeltaBetaWeights.cc:19
DeltaBetaWeights::pfPU_token
edm::EDGetTokenT< edm::View< reco::Candidate > > pfPU_token
Definition: DeltaBetaWeights.h:39
edm::Handle
Definition: AssociativeIterator.h:50
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
DeltaBetaWeights::pfCharged_token
edm::EDGetTokenT< edm::View< reco::Candidate > > pfCharged_token
Definition: DeltaBetaWeights.h:38
pfElectronTranslator_cfi.PFCandidate
PFCandidate
Definition: pfElectronTranslator_cfi.py:6
PVValHelper::eta
Definition: PVValidationHelpers.h:70
DeltaBetaWeights::pfCharged_
edm::InputTag pfCharged_
Definition: DeltaBetaWeights.h:35
DeltaBetaWeights::pfPU_
edm::InputTag pfPU_
Definition: DeltaBetaWeights.h:36
edm::ParameterSet
Definition: ParameterSet.h:47
TrackRefitter_38T_cff.src
src
Definition: TrackRefitter_38T_cff.py:24
cand
Definition: decayParser.h:32
iEvent
int iEvent
Definition: GenABIO.cc:224
reco::LeafCandidate::p4
const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:114
edm::EventSetup
Definition: EventSetup.h:58
reco::Candidate
Definition: Candidate.h:27
DDAxes::phi
DeltaBetaWeights::DeltaBetaWeights
DeltaBetaWeights(const edm::ParameterSet &)
Definition: DeltaBetaWeights.cc:4
reco::LeafCandidate::setP4
void setP4(const LorentzVector &p4) final
set 4-momentum
Definition: LeafCandidate.h:158
eostools.move
def move(src, dest)
Definition: eostools.py:511
HLTMuonOfflineAnalyzer_cfi.deltaR2
deltaR2
Definition: HLTMuonOfflineAnalyzer_cfi.py:105
reco::PFCandidateCollection
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
Definition: PFCandidateFwd.h:12
DeltaBetaWeights::src_
edm::InputTag src_
Definition: DeltaBetaWeights.h:34
reco::PFCandidate::translatePdgIdToType
ParticleType translatePdgIdToType(int pdgid) const
Definition: PFCandidate.cc:209
dqm-mbProfile.log
log
Definition: dqm-mbProfile.py:17
reco::PFCandidate
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
MillePedeFileConverter_cfg.out
out
Definition: MillePedeFileConverter_cfg.py:31
edm::Event
Definition: Event.h:73
DeltaBetaWeights::src_token
edm::EDGetTokenT< edm::View< reco::Candidate > > src_token
Definition: DeltaBetaWeights.h:40
DeltaBetaWeights.h
reco::PFCandidate::setParticleType
void setParticleType(ParticleType type)
set Particle Type
Definition: PFCandidate.cc:256