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