CMS 3D CMS Logo

PuppiProducer.cc
Go to the documentation of this file.
21 
22 #include <memory>
23 
24 // ------------------------------------------------------------------------------------------
26 public:
27  explicit PuppiProducer(const edm::ParameterSet&);
28  ~PuppiProducer() override;
29 
30  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
32  typedef std::vector<LorentzVector> LorentzVectorCollection;
35  typedef std::vector<reco::PFCandidate> PFInputCollection;
36  typedef std::vector<reco::PFCandidate> PFOutputCollection;
37  typedef std::vector<pat::PackedCandidate> PackedOutputCollection;
40 
41 private:
42  virtual void beginJob();
43  void produce(edm::Event&, const edm::EventSetup&) override;
44  virtual void endJob();
45 
72  bool fUseDZ;
74  double fDZCut;
75  double fEtaMinUseDZ;
76  double fPtMaxCharged;
78  double fPtMaxPhotons;
86  double fVtxZCut;
88  std::unique_ptr<PuppiContainer> fPuppiContainer;
89  std::vector<RecoObj> fRecoObjCollection;
90 };
91 
92 // ------------------------------------------------------------------------------------------
94  fPuppiDiagnostics = iConfig.getParameter<bool>("puppiDiagnostics");
95  fPuppiNoLep = iConfig.getParameter<bool>("puppiNoLep");
96  fUseFromPVLooseTight = iConfig.getParameter<bool>("UseFromPVLooseTight");
97  fUseDZ = iConfig.getParameter<bool>("UseDeltaZCut");
98  fUseDZforPileup = iConfig.getParameter<bool>("UseDeltaZCutForPileup");
99  fDZCut = iConfig.getParameter<double>("DeltaZCut");
100  fEtaMinUseDZ = iConfig.getParameter<double>("EtaMinUseDeltaZ");
101  fPtMaxCharged = iConfig.getParameter<double>("PtMaxCharged");
102  fEtaMaxCharged = iConfig.getParameter<double>("EtaMaxCharged");
103  fPtMaxPhotons = iConfig.getParameter<double>("PtMaxPhotons");
104  fEtaMaxPhotons = iConfig.getParameter<double>("EtaMaxPhotons");
105  fNumOfPUVtxsForCharged = iConfig.getParameter<uint>("NumOfPUVtxsForCharged");
106  fDZCutForChargedFromPUVtxs = iConfig.getParameter<double>("DeltaZCutForChargedFromPUVtxs");
107  fUseExistingWeights = iConfig.getParameter<bool>("useExistingWeights");
108  fApplyPhotonProtectionForExistingWeights = iConfig.getParameter<bool>("applyPhotonProtectionForExistingWeights");
109  fClonePackedCands = iConfig.getParameter<bool>("clonePackedCands");
110  fVtxNdofCut = iConfig.getParameter<int>("vtxNdofCut");
111  fVtxZCut = iConfig.getParameter<double>("vtxZCut");
112  fPuppiContainer = std::make_unique<PuppiContainer>(iConfig);
113 
114  tokenPFCandidates_ = consumes<CandidateView>(iConfig.getParameter<edm::InputTag>("candName"));
115  tokenVertices_ = consumes<VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexName"));
116  fUseVertexAssociation = iConfig.getParameter<bool>("useVertexAssociation");
117  vertexAssociationQuality_ = iConfig.getParameter<int>("vertexAssociationQuality");
118  if (fUseVertexAssociation) {
119  tokenVertexAssociation_ = consumes<CandToVertex>(iConfig.getParameter<edm::InputTag>("vertexAssociation"));
121  consumes<edm::ValueMap<int>>(iConfig.getParameter<edm::InputTag>("vertexAssociation"));
122  }
123 
124  fUsePUProxyValue = iConfig.getParameter<bool>("usePUProxyValue");
125 
126  if (fUsePUProxyValue) {
127  puProxyValueToken_ = consumes<double>(iConfig.getParameter<edm::InputTag>("PUProxyValue"));
128  }
129 
130  ptokenPupOut_ = produces<edm::ValueMap<float>>();
131  ptokenP4PupOut_ = produces<edm::ValueMap<LorentzVector>>();
132  ptokenValues_ = produces<edm::ValueMap<reco::CandidatePtr>>();
133 
135  ptokenPackedPuppiCandidates_ = produces<pat::PackedCandidateCollection>();
136  else {
137  ptokenPuppiCandidates_ = produces<reco::PFCandidateCollection>();
138  }
139 
140  if (fPuppiDiagnostics) {
141  ptokenNalgos_ = produces<double>("PuppiNAlgos");
142  ptokenRawAlphas_ = produces<std::vector<double>>("PuppiRawAlphas");
143  ptokenAlphas_ = produces<std::vector<double>>("PuppiAlphas");
144  ptokenAlphasMed_ = produces<std::vector<double>>("PuppiAlphasMed");
145  ptokenAlphasRms_ = produces<std::vector<double>>("PuppiAlphasRms");
146  }
147 }
148 // ------------------------------------------------------------------------------------------
150 // ------------------------------------------------------------------------------------------
152  // Get PFCandidate Collection
153  edm::Handle<CandidateView> hPFProduct;
154  iEvent.getByToken(tokenPFCandidates_, hPFProduct);
155  const CandidateView* pfCol = hPFProduct.product();
156 
157  // Get vertex collection w/PV as the first entry?
159  iEvent.getByToken(tokenVertices_, hVertexProduct);
160  const reco::VertexCollection* pvCol = hVertexProduct.product();
161 
163  edm::ValueMap<int> associationQuality;
165  associatedPV = iEvent.get(tokenVertexAssociation_);
166  associationQuality = iEvent.get(tokenVertexAssociationQuality_);
167  }
168 
169  double puProxyValue = 0.;
170  if (fUsePUProxyValue) {
171  puProxyValue = iEvent.get(puProxyValueToken_);
172  } else {
173  for (auto const& vtx : *pvCol) {
174  if (!vtx.isFake() && vtx.ndof() >= fVtxNdofCut && std::abs(vtx.z()) <= fVtxZCut)
175  ++puProxyValue;
176  }
177  }
178 
179  std::vector<double> lWeights;
180  if (!fUseExistingWeights) {
181  //Fill the reco objects
182  fRecoObjCollection.clear();
183  fRecoObjCollection.reserve(pfCol->size());
184  int iCand = 0;
185  for (auto const& aPF : *pfCol) {
186  RecoObj pReco;
187  pReco.pt = aPF.pt();
188  pReco.eta = aPF.eta();
189  pReco.phi = aPF.phi();
190  pReco.m = aPF.mass();
191  pReco.rapidity = aPF.rapidity();
192  pReco.charge = aPF.charge();
193  pReco.pdgId = aPF.pdgId();
194  const reco::Vertex* closestVtx = nullptr;
195  double pDZ = -9999;
196  double pD0 = -9999;
197  uint pVtxId = 0;
198  bool isLepton = ((std::abs(pReco.pdgId) == 11) || (std::abs(pReco.pdgId) == 13));
199  const pat::PackedCandidate* lPack = dynamic_cast<const pat::PackedCandidate*>(&aPF);
200 
201  if (fUseVertexAssociation) {
202  const reco::VertexRef& PVOrig = associatedPV[reco::CandidatePtr(hPFProduct, iCand)];
203  int quality = associationQuality[reco::CandidatePtr(hPFProduct, iCand)];
204  if (PVOrig.isNonnull() && (quality >= vertexAssociationQuality_)) {
205  closestVtx = PVOrig.get();
206  pVtxId = PVOrig.key();
207  }
208  if (std::abs(pReco.charge) == 0)
209  pReco.id = 0;
210  else if (fPuppiNoLep && isLepton)
211  pReco.id = 3;
212  else if (closestVtx != nullptr && pVtxId == 0)
213  pReco.id = 1; // Associated to main vertex
214  else if (closestVtx != nullptr && pVtxId > 0)
215  pReco.id = 2; // Associated to PU
216  else
217  pReco.id = 0; // Unassociated
218  } else if (lPack == nullptr) {
219  const reco::PFCandidate* pPF = dynamic_cast<const reco::PFCandidate*>(&aPF);
220  double curdz = 9999;
221  int closestVtxForUnassociateds = -9999;
222  const reco::TrackRef aTrackRef = pPF->trackRef();
223  bool lFirst = true;
224  for (auto const& aV : *pvCol) {
225  if (lFirst) {
226  if (aTrackRef.isNonnull()) {
227  pDZ = aTrackRef->dz(aV.position());
228  pD0 = aTrackRef->d0();
229  } else if (pPF->gsfTrackRef().isNonnull()) {
230  pDZ = pPF->gsfTrackRef()->dz(aV.position());
231  pD0 = pPF->gsfTrackRef()->d0();
232  }
233  lFirst = false;
234  if (pDZ > -9999)
235  pVtxId = 0;
236  }
237  if (aTrackRef.isNonnull() && aV.trackWeight(pPF->trackRef()) > 0) {
238  closestVtx = &aV;
239  break;
240  }
241  // in case it's unassocciated, keep more info
242  double tmpdz = 99999;
243  if (aTrackRef.isNonnull())
244  tmpdz = aTrackRef->dz(aV.position());
245  else if (pPF->gsfTrackRef().isNonnull())
246  tmpdz = pPF->gsfTrackRef()->dz(aV.position());
247  if (std::abs(tmpdz) < curdz) {
248  curdz = std::abs(tmpdz);
249  closestVtxForUnassociateds = pVtxId;
250  }
251  pVtxId++;
252  }
253  int tmpFromPV = 0;
254  // mocking the miniAOD definitions
255  if (std::abs(pReco.charge) > 0) {
256  if (closestVtx != nullptr && pVtxId > 0)
257  tmpFromPV = 0;
258  if (closestVtx != nullptr && pVtxId == 0)
259  tmpFromPV = 3;
260  if (closestVtx == nullptr && closestVtxForUnassociateds == 0)
261  tmpFromPV = 2;
262  if (closestVtx == nullptr && closestVtxForUnassociateds != 0)
263  tmpFromPV = 1;
264  }
265  pReco.dZ = pDZ;
266  pReco.d0 = pD0;
267  pReco.id = 0;
268  if (std::abs(pReco.charge) == 0) {
269  pReco.id = 0;
270  } else {
271  if (fPuppiNoLep && isLepton)
272  pReco.id = 3;
273  else if (tmpFromPV == 0) {
274  pReco.id = 2;
275  if (fNumOfPUVtxsForCharged > 0 and (pVtxId <= fNumOfPUVtxsForCharged) and
277  pReco.id = 1;
278  } else if (tmpFromPV == 3)
279  pReco.id = 1;
280  else if (tmpFromPV == 1 || tmpFromPV == 2) {
281  pReco.id = 0;
282  if ((fPtMaxCharged > 0) and (pReco.pt > fPtMaxCharged))
283  pReco.id = 1;
284  else if (std::abs(pReco.eta) > fEtaMaxCharged)
285  pReco.id = 1;
286  else if ((fUseDZ) && (std::abs(pReco.eta) >= fEtaMinUseDZ) && (std::abs(pDZ) < fDZCut))
287  pReco.id = 1;
288  else if ((fUseDZforPileup) && (std::abs(pReco.eta) >= fEtaMinUseDZ) && (std::abs(pDZ) >= fDZCut))
289  pReco.id = 2;
290  else if (fUseFromPVLooseTight && tmpFromPV == 1)
291  pReco.id = 2;
292  else if (fUseFromPVLooseTight && tmpFromPV == 2)
293  pReco.id = 1;
294  }
295  }
296  } else if (lPack->vertexRef().isNonnull()) {
297  pDZ = lPack->dz();
298  pD0 = lPack->dxy();
299  pReco.dZ = pDZ;
300  pReco.d0 = pD0;
301 
302  pReco.id = 0;
303  if (std::abs(pReco.charge) == 0) {
304  pReco.id = 0;
305  }
306  if (std::abs(pReco.charge) > 0) {
307  if (fPuppiNoLep && isLepton) {
308  pReco.id = 3;
309  } else if (lPack->fromPV() == 0) {
310  pReco.id = 2;
312  for (size_t puVtx_idx = 1; puVtx_idx <= fNumOfPUVtxsForCharged && puVtx_idx < pvCol->size();
313  ++puVtx_idx) {
314  if (lPack->fromPV(puVtx_idx) >= 2) {
315  pReco.id = 1;
316  break;
317  }
318  }
319  }
320  } else if (lPack->fromPV() == (pat::PackedCandidate::PVUsedInFit)) {
321  pReco.id = 1;
322  } else if (lPack->fromPV() == (pat::PackedCandidate::PVTight) ||
323  lPack->fromPV() == (pat::PackedCandidate::PVLoose)) {
324  pReco.id = 0;
325  if ((fPtMaxCharged > 0) and (pReco.pt > fPtMaxCharged))
326  pReco.id = 1;
327  else if (std::abs(pReco.eta) > fEtaMaxCharged)
328  pReco.id = 1;
329  else if ((fUseDZ) && (std::abs(pReco.eta) >= fEtaMinUseDZ) && (std::abs(pDZ) < fDZCut))
330  pReco.id = 1;
331  else if ((fUseDZforPileup) && (std::abs(pReco.eta) >= fEtaMinUseDZ) && (std::abs(pDZ) >= fDZCut))
332  pReco.id = 2;
334  pReco.id = 2;
336  pReco.id = 1;
337  }
338  }
339  }
340 
341  fRecoObjCollection.push_back(pReco);
342  iCand++;
343  }
344 
345  fPuppiContainer->initialize(fRecoObjCollection);
346  fPuppiContainer->setPUProxy(puProxyValue);
347 
348  //Compute the weights and get the particles
349  lWeights = fPuppiContainer->puppiWeights();
350  } else {
351  //Use the existing weights
352  int lPackCtr = 0;
353  lWeights.reserve(pfCol->size());
354  for (auto const& aPF : *pfCol) {
355  const pat::PackedCandidate* lPack = dynamic_cast<const pat::PackedCandidate*>(&aPF);
356  float curpupweight = -1.;
357  if (lPack == nullptr) {
358  // throw error
360  "PuppiProducer: cannot get weights since inputs are not PackedCandidates");
361  } else {
362  if (fPuppiNoLep) {
363  curpupweight = lPack->puppiWeightNoLep();
364  } else {
365  curpupweight = lPack->puppiWeight();
366  }
367  }
368  // Optional: Protect high pT photons (important for gamma to hadronic recoil balance) for existing weights.
369  if (fApplyPhotonProtectionForExistingWeights && (fPtMaxPhotons > 0) && (lPack->pdgId() == 22) &&
370  (std::abs(lPack->eta()) < fEtaMaxPhotons) && (lPack->pt() > fPtMaxPhotons))
371  curpupweight = 1;
372  lWeights.push_back(curpupweight);
373  lPackCtr++;
374  }
375  }
376 
377  //Fill it into the event
378  edm::ValueMap<float> lPupOut;
379  edm::ValueMap<float>::Filler lPupFiller(lPupOut);
380  lPupFiller.insert(hPFProduct, lWeights.begin(), lWeights.end());
381  lPupFiller.fill();
382 
383  // This is a dummy to access the "translate" method which is a
384  // non-static member function even though it doesn't need to be.
385  // Will fix in the future.
386  static const reco::PFCandidate dummySinceTranslateIsNotStatic;
387 
388  // Fill a new PF/Packed Candidate Collection and write out the ValueMap of the new p4s.
389  // Since the size of the ValueMap must be equal to the input collection, we need
390  // to search the "puppi" particles to find a match for each input. If none is found,
391  // the input is set to have a four-vector of 0,0,0,0
392  PFOutputCollection fPuppiCandidates;
393  PackedOutputCollection fPackedPuppiCandidates;
394 
396  LorentzVectorCollection puppiP4s;
397  std::vector<reco::CandidatePtr> values(hPFProduct->size());
398 
399  int iCand = -1;
400  puppiP4s.reserve(hPFProduct->size());
402  fPackedPuppiCandidates.reserve(hPFProduct->size());
403  else
404  fPuppiCandidates.reserve(hPFProduct->size());
405  for (auto const& aCand : *hPFProduct) {
406  ++iCand;
407  std::unique_ptr<pat::PackedCandidate> pCand;
408  std::unique_ptr<reco::PFCandidate> pfCand;
409 
411  const pat::PackedCandidate* cand = dynamic_cast<const pat::PackedCandidate*>(&aCand);
412  if (!cand)
413  throw edm::Exception(edm::errors::LogicError, "PuppiProducer: inputs are not PackedCandidates");
414  pCand = std::make_unique<pat::PackedCandidate>(*cand);
415  } else {
416  auto id = dummySinceTranslateIsNotStatic.translatePdgIdToType(aCand.pdgId());
417  const reco::PFCandidate* cand = dynamic_cast<const reco::PFCandidate*>(&aCand);
418  pfCand = std::make_unique<reco::PFCandidate>(cand ? *cand : reco::PFCandidate(aCand.charge(), aCand.p4(), id));
419  }
420 
421  // Here, we are using new weights computed and putting them in the packed candidates.
423  if (fPuppiNoLep)
424  pCand->setPuppiWeight(pCand->puppiWeight(), lWeights[iCand]);
425  else
426  pCand->setPuppiWeight(lWeights[iCand], pCand->puppiWeightNoLep());
427  }
428 
429  puppiP4s.emplace_back(lWeights[iCand] * aCand.px(),
430  lWeights[iCand] * aCand.py(),
431  lWeights[iCand] * aCand.pz(),
432  lWeights[iCand] * aCand.energy());
433 
434  // Here, we are using existing weights, or we're using packed candidates.
435  // That is, whether or not we recomputed the weights, we store the
436  // source candidate appropriately, and set the p4 of the packed candidate.
438  pCand->setP4(puppiP4s.back());
439  pCand->setSourceCandidatePtr(aCand.sourceCandidatePtr(0));
440  fPackedPuppiCandidates.push_back(*pCand);
441  } else {
442  pfCand->setP4(puppiP4s.back());
443  pfCand->setSourceCandidatePtr(aCand.sourceCandidatePtr(0));
444  fPuppiCandidates.push_back(*pfCand);
445  }
446  }
447 
448  //Compute the modified p4s
449  edm::ValueMap<LorentzVector>::Filler p4PupFiller(p4PupOut);
450  p4PupFiller.insert(hPFProduct, puppiP4s.begin(), puppiP4s.end());
451  p4PupFiller.fill();
452 
453  iEvent.emplace(ptokenPupOut_, lPupOut);
454  iEvent.emplace(ptokenP4PupOut_, p4PupOut);
457  iEvent.emplace(ptokenPackedPuppiCandidates_, fPackedPuppiCandidates);
458  for (unsigned int ic = 0, nc = oh->size(); ic < nc; ++ic) {
459  reco::CandidatePtr pkref(oh, ic);
460  values[ic] = pkref;
461  }
462  } else {
464  for (unsigned int ic = 0, nc = oh->size(); ic < nc; ++ic) {
465  reco::CandidatePtr pkref(oh, ic);
466  values[ic] = pkref;
467  }
468  }
471  filler.insert(hPFProduct, values.begin(), values.end());
472  filler.fill();
473  iEvent.emplace(ptokenValues_, pfMap_p);
474 
477  // all the different alphas per particle
478  // THE alpha per particle
479  std::vector<double> theAlphas(fPuppiContainer->puppiAlphas());
480  std::vector<double> theAlphasMed(fPuppiContainer->puppiAlphasMed());
481  std::vector<double> theAlphasRms(fPuppiContainer->puppiAlphasRMS());
482  std::vector<double> alphas(fPuppiContainer->puppiRawAlphas());
483  double nalgos(fPuppiContainer->puppiNAlgos());
484 
485  iEvent.emplace(ptokenRawAlphas_, alphas);
486  iEvent.emplace(ptokenNalgos_, nalgos);
487  iEvent.emplace(ptokenAlphas_, theAlphas);
488  iEvent.emplace(ptokenAlphasMed_, theAlphasMed);
489  iEvent.emplace(ptokenAlphasRms_, theAlphasRms);
490  }
491 }
492 
493 // ------------------------------------------------------------------------------------------
495 // ------------------------------------------------------------------------------------------
497 // ------------------------------------------------------------------------------------------
500  desc.add<bool>("puppiDiagnostics", false);
501  desc.add<bool>("puppiNoLep", false);
502  desc.add<bool>("UseFromPVLooseTight", false);
503  desc.add<bool>("UseDeltaZCut", true);
504  desc.add<bool>("UseDeltaZCutForPileup", true);
505  desc.add<double>("DeltaZCut", 0.3);
506  desc.add<double>("EtaMinUseDeltaZ", 0.);
507  desc.add<double>("PtMaxCharged", -1.);
508  desc.add<double>("EtaMaxCharged", 99999.);
509  desc.add<double>("PtMaxPhotons", -1.);
510  desc.add<double>("EtaMaxPhotons", 2.5);
511  desc.add<double>("PtMaxNeutrals", 200.);
512  desc.add<double>("PtMaxNeutralsStartSlope", 0.);
513  desc.add<uint>("NumOfPUVtxsForCharged", 0);
514  desc.add<double>("DeltaZCutForChargedFromPUVtxs", 0.2);
515  desc.add<bool>("useExistingWeights", false);
516  desc.add<bool>("applyPhotonProtectionForExistingWeights", false);
517  desc.add<bool>("clonePackedCands", false);
518  desc.add<int>("vtxNdofCut", 4);
519  desc.add<double>("vtxZCut", 24);
520  desc.add<edm::InputTag>("candName", edm::InputTag("particleFlow"));
521  desc.add<edm::InputTag>("vertexName", edm::InputTag("offlinePrimaryVertices"));
522  desc.add<bool>("useVertexAssociation", false);
523  desc.add<int>("vertexAssociationQuality", 0);
524  desc.add<edm::InputTag>("vertexAssociation", edm::InputTag(""));
525  desc.add<bool>("applyCHS", true);
526  desc.add<bool>("invertPuppi", false);
527  desc.add<bool>("useExp", false);
528  desc.add<double>("MinPuppiWeight", .01);
529  desc.add<bool>("usePUProxyValue", false);
530  desc.add<edm::InputTag>("PUProxyValue", edm::InputTag(""));
531 
533 
534  descriptions.add("PuppiProducer", desc);
535 }
536 //define this as a plug-in
float puppiWeight() const
std::vector< reco::PFCandidate > PFInputCollection
reco::GsfTrackRef gsfTrackRef() const
Definition: PFCandidate.cc:467
bool fUseFromPVLooseTight
float puppiWeightNoLep() const
Weight from full PUPPI.
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDPutTokenT< pat::PackedCandidateCollection > ptokenPackedPuppiCandidates_
std::string fPFName
edm::EDPutTokenT< std::vector< double > > ptokenAlphasMed_
edm::EDGetTokenT< PuppiContainer > tokenPuppiContainer_
float d0
Definition: RecoObj.h:37
double fDZCutForChargedFromPUVtxs
virtual void endJob()
static void fillDescriptionsPuppiAlgo(edm::ParameterSetDescription &desc)
Definition: PuppiAlgo.cc:218
Definition: RecoObj.h:4
bool fPuppiDiagnostics
T const * product() const
Definition: Handle.h:70
bool fUseVertexAssociation
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:53
std::vector< RecoObj > fRecoObjCollection
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
std::vector< LorentzVector > LorentzVectorCollection
std::string fPVName
bool isLepton(const Candidate &part)
Definition: pdgIdUtils.h:13
const PVAssoc fromPV(size_t ipv=0) const
edm::EDGetTokenT< CandToVertex > tokenVertexAssociation_
virtual void beginJob()
edm::EDPutTokenT< reco::PFCandidateCollection > ptokenPuppiCandidates_
edm::EDPutTokenT< double > ptokenNalgos_
int pdgId() const override
PDG identifier.
math::XYZTLorentzVector LorentzVector
const reco::VertexRef vertexRef() const
double fEtaMaxCharged
std::unique_ptr< PuppiContainer > fPuppiContainer
int charge
Definition: RecoObj.h:38
key_type key() const
Accessor for product key.
Definition: Ref.h:250
edm::EDGetTokenT< PFOutputCollection > tokenPuppiCandidates_
float dZ
Definition: RecoObj.h:36
size_type size() const
edm::EDPutTokenT< edm::ValueMap< reco::CandidatePtr > > ptokenValues_
float rapidity
Definition: RecoObj.h:26
int vertexAssociationQuality_
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
reco::VertexCollection VertexCollection
string quality
edm::EDGetTokenT< double > puProxyValueToken_
float pt
Definition: RecoObj.h:26
int iEvent
Definition: GenABIO.cc:224
double eta() const override
momentum pseudorapidity
int id
Definition: RecoObj.h:27
float phi
Definition: RecoObj.h:26
std::vector< reco::PFCandidate > PFOutputCollection
edm::EDGetTokenT< edm::ValueMap< int > > tokenVertexAssociationQuality_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool fClonePackedCands
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::View< reco::Candidate > CandidateView
void produce(edm::Event &, const edm::EventSetup &) override
std::string fPuppiName
std::vector< pat::PackedCandidate > PackedOutputCollection
PuppiProducer(const edm::ParameterSet &)
edm::EDGetTokenT< VertexCollection > tokenVertices_
edm::EDPutTokenT< edm::ValueMap< float > > ptokenPupOut_
edm::EDPutTokenT< std::vector< double > > ptokenAlphas_
edm::EDPutTokenT< std::vector< double > > ptokenRawAlphas_
edm::EDGetTokenT< CandidateView > tokenPFCandidates_
edm::Association< reco::VertexCollection > CandToVertex
edm::Ptr< Candidate > CandidatePtr
persistent reference to an object in a collection of Candidate objects
Definition: CandidateFwd.h:25
float m
Definition: RecoObj.h:26
double pt() const override
transverse momentum
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
double fPtMaxPhotons
double fPtMaxCharged
void add(std::string const &label, ParameterSetDescription const &psetDescription)
int pdgId
Definition: RecoObj.h:28
float eta
Definition: RecoObj.h:26
bool fApplyPhotonProtectionForExistingWeights
uint fNumOfPUVtxsForCharged
edm::EDPutTokenT< edm::ValueMap< LorentzVector > > ptokenP4PupOut_
bool fUseExistingWeights
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
ParticleType translatePdgIdToType(int pdgid) const
Definition: PFCandidate.cc:231
edm::View< reco::PFCandidate > PFView
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:232
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:435
virtual float dxy() const
dxy with respect to the PV ref
edm::EDPutTokenT< std::vector< double > > ptokenAlphasRms_
double fEtaMaxPhotons
edm::EDGetTokenT< PackedOutputCollection > tokenPackedPuppiCandidates_
virtual float dz(size_t ipv=0) const
dz with respect to the PV[ipv]
~PuppiProducer() override
double fEtaMinUseDZ