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