CMS 3D CMS Logo

PuppiPhoton.cc
Go to the documentation of this file.
22 
23 #include <memory>
24 
25 // ------------------------------------------------------------------------------------------
27 public:
28  explicit PuppiPhoton(const edm::ParameterSet &);
29  ~PuppiPhoton() override;
30 
31  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
33  typedef std::vector<LorentzVector> LorentzVectorCollection;
35  typedef std::vector<reco::PFCandidate> PFOutputCollection;
37 
38 private:
39  void produce(edm::Event &, const edm::EventSetup &) override;
40  bool matchPFCandidate(const reco::Candidate *iPF, const reco::Candidate *iPho);
47  double pt_;
48  double eta_;
49  bool usePFRef_;
53  std::vector<double> dRMatch_;
54  std::vector<int32_t> pdgIds_;
55  std::unique_ptr<PFOutputCollection> corrCandidates_;
56  double weight_;
58 };
59 
60 // ------------------------------------------------------------------------------------------
62  tokenPFCandidates_ = consumes<CandidateView>(iConfig.getParameter<edm::InputTag>("candName"));
63  tokenPuppiCandidates_ = consumes<CandidateView>(iConfig.getParameter<edm::InputTag>("puppiCandName"));
64  usePFphotons_ = iConfig.getParameter<bool>("usePFphotons");
65  if (!usePFphotons_)
66  tokenPhotonCandidates_ = consumes<CandidateView>(iConfig.getParameter<edm::InputTag>("photonName"));
67  usePhotonId_ = !(iConfig.getParameter<edm::InputTag>("photonId")).label().empty();
68  if (usePhotonId_)
69  tokenPhotonId_ = consumes<edm::ValueMap<bool>>(iConfig.getParameter<edm::InputTag>("photonId"));
70  runOnMiniAOD_ = iConfig.getParameter<bool>("runOnMiniAOD");
71  if (!runOnMiniAOD_)
72  reco2pf_ =
73  consumes<edm::ValueMap<std::vector<reco::PFCandidateRef>>>(iConfig.getParameter<edm::InputTag>("recoToPFMap"));
74  useValueMap_ = iConfig.getParameter<bool>("useValueMap");
75  if (useValueMap_)
76  tokenWeights_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("weightsName"));
77 
78  pt_ = iConfig.getParameter<double>("pt");
79  eta_ = iConfig.getParameter<double>("eta");
80  dRMatch_ = iConfig.getParameter<std::vector<double>>("dRMatch");
81  pdgIds_ = iConfig.getParameter<std::vector<int32_t>>("pdgids");
82  usePFRef_ = iConfig.getParameter<bool>("useRefs");
83  weight_ = iConfig.getParameter<double>("weight");
84  produces<PFOutputCollection>();
85  produces<edm::ValueMap<reco::CandidatePtr>>();
86 }
87 // ------------------------------------------------------------------------------------------
89 // ------------------------------------------------------------------------------------------
91  int iC = -1;
92  std::vector<const reco::Candidate *> phoCands;
93  std::vector<uint16_t> phoIndx;
94 
96  if (!runOnMiniAOD_)
97  iEvent.getByToken(reco2pf_, reco2pf);
98 
99  // Get PFCandidate Collection
100  edm::Handle<CandidateView> hPFProduct;
101  iEvent.getByToken(tokenPFCandidates_, hPFProduct);
102  const CandidateView *pfCol = hPFProduct.product();
103 
104  edm::Handle<CandidateView> hPuppiProduct;
105  iEvent.getByToken(tokenPuppiCandidates_, hPuppiProduct);
106  const CandidateView *pupCol = hPuppiProduct.product();
107  if (usePFphotons_) {
108  for (const auto &pho : *pfCol) {
109  iC++;
110  if (pho.pt() < pt_)
111  continue;
112  if (std::abs(pho.pdgId()) != 22)
113  continue;
114  if (fabs(pho.eta()) < eta_) {
115  phoIndx.push_back(iC);
116  phoCands.push_back(&pho);
117  }
118  }
119  } else {
120  edm::Handle<CandidateView> hPhoProduct;
121  iEvent.getByToken(tokenPhotonCandidates_, hPhoProduct);
122  const CandidateView *phoCol = hPhoProduct.product();
123 
125  if (usePhotonId_)
126  iEvent.getByToken(tokenPhotonId_, photonId);
127 
128  for (CandidateView::const_iterator itPho = phoCol->begin(); itPho != phoCol->end(); itPho++) {
129  iC++;
130  bool passObject = false;
131  if (itPho->isPhoton() && usePhotonId_)
132  passObject = (*photonId)[phoCol->ptrAt(iC)];
133  if (itPho->pt() < pt_)
134  continue;
135  if (!passObject && usePhotonId_)
136  continue;
137  if (runOnMiniAOD_) {
138  const pat::Photon *pPho = dynamic_cast<const pat::Photon *>(&(*itPho));
139  if (pPho != nullptr) {
141  if (fabs(ref->eta()) < eta_) {
142  phoIndx.push_back(ref.key());
143  phoCands.push_back(&(*(pfCol->ptrAt(ref.key()))));
144  }
145  }
146  continue;
147  }
148  const pat::Electron *pElectron = dynamic_cast<const pat::Electron *>(&(*itPho));
149  if (pElectron != nullptr) {
151  if (fabs(ref->eta()) < eta_) {
152  phoIndx.push_back(ref.key());
153  phoCands.push_back(&(*(pfCol->ptrAt(ref.key()))));
154  }
155  }
156  } else {
157  for (const edm::Ref<std::vector<reco::PFCandidate>> &ref : (*reco2pf)[phoCol->ptrAt(iC)]) {
158  if (fabs(ref->eta()) < eta_) {
159  phoIndx.push_back(ref.key());
160  phoCands.push_back(&(*(pfCol->ptrAt(ref.key()))));
161  }
162  }
163  }
164  }
165  }
166  //Get Weights
168  if (useValueMap_)
169  iEvent.getByToken(tokenWeights_, pupWeights);
170  std::unique_ptr<edm::ValueMap<LorentzVector>> p4PupOut(new edm::ValueMap<LorentzVector>());
171  LorentzVectorCollection puppiP4s;
172  std::vector<reco::CandidatePtr> values(hPFProduct->size());
173  int iPF = 0;
174  std::vector<float> lWeights;
175  static const reco::PFCandidate dummySinceTranslateIsNotStatic;
176  corrCandidates_ = std::make_unique<PFOutputCollection>();
177  std::set<int> foundPhoIndex;
178  for (CandidateView::const_iterator itPF = pupCol->begin(); itPF != pupCol->end(); itPF++) {
179  auto id = dummySinceTranslateIsNotStatic.translatePdgIdToType(itPF->pdgId());
180  const reco::PFCandidate *pPF = dynamic_cast<const reco::PFCandidate *>(&(*itPF));
181  reco::PFCandidate pCand(pPF ? *pPF : reco::PFCandidate(itPF->charge(), itPF->p4(), id));
182  LorentzVector pVec = itPF->p4();
183  float pWeight = 1.;
184  if (useValueMap_)
185  pWeight = (*pupWeights)[pupCol->ptrAt(iPF)];
186  if (!usePFRef_) {
187  int iPho = -1;
188  for (std::vector<const reco::Candidate *>::iterator itPho = phoCands.begin(); itPho != phoCands.end(); itPho++) {
189  iPho++;
190  if ((!matchPFCandidate(&(*itPF), *itPho)) || (foundPhoIndex.count(iPho) != 0))
191  continue;
192  pWeight = weight_;
193  if (!useValueMap_ && itPF->pt() != 0)
194  pWeight = pWeight * (phoCands[iPho]->pt() / itPF->pt());
195  if (!useValueMap_ && itPF->pt() == 0)
196  pVec.SetPxPyPzE(phoCands[iPho]->px() * pWeight,
197  phoCands[iPho]->py() * pWeight,
198  phoCands[iPho]->pz() * pWeight,
199  phoCands[iPho]->energy() * pWeight);
200  foundPhoIndex.insert(iPho);
201  }
202  } else {
203  int iPho = -1;
204  for (std::vector<uint16_t>::const_iterator itPho = phoIndx.begin(); itPho != phoIndx.end(); itPho++) {
205  iPho++;
206  if (pupCol->refAt(iPF).key() != *itPho)
207  continue;
208  pWeight = weight_;
209  if (!useValueMap_ && itPF->pt() != 0)
210  pWeight = pWeight * (phoCands[iPho]->pt() / itPF->pt());
211  if (!useValueMap_ && itPF->pt() == 0)
212  pVec.SetPxPyPzE(phoCands[iPho]->px() * pWeight,
213  phoCands[iPho]->py() * pWeight,
214  phoCands[iPho]->pz() * pWeight,
215  phoCands[iPho]->energy() * pWeight);
216  foundPhoIndex.insert(iPho);
217  }
218  }
219  if (itPF->pt() != 0)
220  pVec.SetPxPyPzE(itPF->px() * pWeight, itPF->py() * pWeight, itPF->pz() * pWeight, itPF->energy() * pWeight);
221 
222  lWeights.push_back(pWeight);
223  pCand.setP4(pVec);
224  puppiP4s.push_back(pVec);
225  pCand.setSourceCandidatePtr(itPF->sourceCandidatePtr(0));
226  corrCandidates_->push_back(pCand);
227  iPF++;
228  }
229  //Add the missing pfcandidates
230  for (unsigned int iPho = 0; iPho < phoCands.size(); iPho++) {
231  if (foundPhoIndex.count(iPho) != 0)
232  continue;
233  auto id = dummySinceTranslateIsNotStatic.translatePdgIdToType(phoCands[iPho]->pdgId());
234  reco::PFCandidate pCand(reco::PFCandidate(phoCands[iPho]->charge(), phoCands[iPho]->p4(), id));
235  pCand.setSourceCandidatePtr(phoCands[iPho]->sourceCandidatePtr(0));
236  LorentzVector pVec = phoCands[iPho]->p4();
237  pVec.SetPxPyPzE(phoCands[iPho]->px() * weight_,
238  phoCands[iPho]->py() * weight_,
239  phoCands[iPho]->pz() * weight_,
240  phoCands[iPho]->energy() * weight_);
241  pCand.setP4(pVec);
242  lWeights.push_back(weight_);
243  puppiP4s.push_back(pVec);
244  corrCandidates_->push_back(pCand);
245  }
246  //Fill it into the event
248  for (unsigned int ic = 0, nc = pupCol->size(); ic < nc; ++ic) {
249  reco::CandidatePtr pkref(oh, ic);
250  values[ic] = pkref;
251  }
252  std::unique_ptr<edm::ValueMap<reco::CandidatePtr>> pfMap_p(new edm::ValueMap<reco::CandidatePtr>());
254  filler.insert(hPFProduct, values.begin(), values.end());
255  filler.fill();
256  iEvent.put(std::move(pfMap_p));
257 }
258 // ------------------------------------------------------------------------------------------
260  if (iPF->pdgId() != iPho->pdgId())
261  return false;
262  double lDR = deltaR(iPF->eta(), iPF->phi(), iPho->eta(), iPho->phi());
263  for (unsigned int i0 = 0; i0 < pdgIds_.size(); i0++) {
264  if (std::abs(iPF->pdgId()) == pdgIds_[i0] && lDR < dRMatch_[i0])
265  return true;
266  }
267  return false;
268 }
269 // ------------------------------------------------------------------------------------------
271  //The following says we do not know what parameters are allowed so do no validation
272  // Please change this to state exactly what you do use, even if it is no parameters
274  desc.setUnknown();
275  descriptions.addDefault(desc);
276 }
277 //define this as a plug-in
std::vector< int32_t > pdgIds_
Definition: PuppiPhoton.cc:54
double pt_
Definition: PuppiPhoton.cc:47
Analysis-level Photon class.
Definition: Photon.h:46
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
PuppiPhoton(const edm::ParameterSet &)
Definition: PuppiPhoton.cc:61
RefToBase< value_type > refAt(size_type i) const
Ptr< value_type > ptrAt(size_type i) const
~PuppiPhoton() override
Definition: PuppiPhoton.cc:88
bool runOnMiniAOD_
Definition: PuppiPhoton.cc:51
T const * product() const
Definition: Handle.h:70
edm::View< reco::PFCandidate > PFView
Definition: PuppiPhoton.cc:36
edm::View< reco::Candidate > CandidateView
Definition: PuppiPhoton.cc:34
void setSourceCandidatePtr(const PFCandidatePtr &ptr)
Definition: PFCandidate.h:120
bool usePFRef_
Definition: PuppiPhoton.cc:49
edm::RefVector< pat::PackedCandidateCollection > associatedPackedPFCandidates() const
References to PFCandidates linked to this object (e.g. for isolation vetos or masking before jet recl...
size_type size() const
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
bool usePFphotons_
Definition: PuppiPhoton.cc:50
int iEvent
Definition: GenABIO.cc:224
void addDefault(ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< edm::ValueMap< float > > tokenWeights_
Definition: PuppiPhoton.cc:45
edm::EDGetTokenT< edm::ValueMap< bool > > tokenPhotonId_
Definition: PuppiPhoton.cc:46
edm::EDGetTokenT< CandidateView > tokenPFCandidates_
Definition: PuppiPhoton.cc:41
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
bool matchPFCandidate(const reco::Candidate *iPF, const reco::Candidate *iPho)
Definition: PuppiPhoton.cc:259
bool usePhotonId_
Definition: PuppiPhoton.cc:52
double weight_
Definition: PuppiPhoton.cc:56
virtual int pdgId() const =0
PDG identifier.
double eta_
Definition: PuppiPhoton.cc:48
std::vector< reco::PFCandidate > PFOutputCollection
Definition: PuppiPhoton.cc:35
void produce(edm::Event &, const edm::EventSetup &) override
Definition: PuppiPhoton.cc:90
Analysis-level electron class.
Definition: Electron.h:51
edm::EDGetTokenT< CandidateView > tokenPuppiCandidates_
Definition: PuppiPhoton.cc:42
std::unique_ptr< PFOutputCollection > corrCandidates_
Definition: PuppiPhoton.cc:55
std::vector< double > dRMatch_
Definition: PuppiPhoton.cc:53
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
ParticleType translatePdgIdToType(int pdgid) const
Definition: PFCandidate.cc:231
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:88
const_iterator begin() const
math::XYZTLorentzVector LorentzVector
Definition: PuppiPhoton.cc:32
edm::EDGetTokenT< edm::ValueMap< std::vector< reco::PFCandidateRef > > > reco2pf_
Definition: PuppiPhoton.cc:44
std::vector< LorentzVector > LorentzVectorCollection
Definition: PuppiPhoton.cc:33
const_iterator end() const
void setP4(const LorentzVector &p4) final
set 4-momentum
def move(src, dest)
Definition: eostools.py:511
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: PuppiPhoton.cc:270
bool useValueMap_
Definition: PuppiPhoton.cc:57
edm::RefVector< pat::PackedCandidateCollection > associatedPackedPFCandidates() const
References to PFCandidates linked to this object (e.g. for isolation vetos or masking before jet recl...
virtual double phi() const =0
momentum azimuthal angle
virtual double eta() const =0
momentum pseudorapidity
edm::EDGetTokenT< CandidateView > tokenPhotonCandidates_
Definition: PuppiPhoton.cc:43