CMS 3D CMS Logo

PuppiProducer.cc
Go to the documentation of this file.
1 // system include files
2 #include <memory>
3 
4 // user include files
19 //Main File
22 
23 
24 // ------------------------------------------------------------------------------------------
26  fPuppiDiagnostics = iConfig.getParameter<bool>("puppiDiagnostics");
27  fPuppiForLeptons = iConfig.getParameter<bool>("puppiForLeptons");
28  fUseFromPVLooseTight = iConfig.getParameter<bool>("UseFromPVLooseTight");
29  fUseDZ = iConfig.getParameter<bool>("UseDeltaZCut");
30  fDZCut = iConfig.getParameter<double>("DeltaZCut");
31  fPtMaxCharged = iConfig.getParameter<double>("PtMaxCharged");
32  fEtaMaxCharged = iConfig.getParameter<double>("EtaMaxCharged");
33  fUseExistingWeights = iConfig.getParameter<bool>("useExistingWeights");
34  fUseWeightsNoLep = iConfig.getParameter<bool>("useWeightsNoLep");
35  fClonePackedCands = iConfig.getParameter<bool>("clonePackedCands");
36  fVtxNdofCut = iConfig.getParameter<int>("vtxNdofCut");
37  fVtxZCut = iConfig.getParameter<double>("vtxZCut");
38  fPuppiContainer = std::unique_ptr<PuppiContainer> ( new PuppiContainer(iConfig) );
39 
41  = consumes<CandidateView>(iConfig.getParameter<edm::InputTag>("candName"));
43  = consumes<VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexName"));
44 
45 
46  produces<edm::ValueMap<float> > ();
47  produces<edm::ValueMap<LorentzVector> > ();
48  produces< edm::ValueMap<reco::CandidatePtr> >();
49 
50  if (fUseExistingWeights || fClonePackedCands)
51  produces<pat::PackedCandidateCollection>();
52  else
53  produces<reco::PFCandidateCollection>();
54 
55  if (fPuppiDiagnostics){
56  produces<double> ("PuppiNAlgos");
57  produces<std::vector<double>> ("PuppiRawAlphas");
58  produces<std::vector<double>> ("PuppiAlphas");
59  produces<std::vector<double>> ("PuppiAlphasMed");
60  produces<std::vector<double>> ("PuppiAlphasRms");
61  }
62 }
63 // ------------------------------------------------------------------------------------------
65 }
66 // ------------------------------------------------------------------------------------------
68 
69  // Get PFCandidate Collection
70  edm::Handle<CandidateView> hPFProduct;
71  iEvent.getByToken(tokenPFCandidates_,hPFProduct);
72  const CandidateView *pfCol = hPFProduct.product();
73 
74  // Get vertex collection w/PV as the first entry?
76  iEvent.getByToken(tokenVertices_,hVertexProduct);
77  const reco::VertexCollection *pvCol = hVertexProduct.product();
78 
79  int npv = 0;
80  const reco::VertexCollection::const_iterator vtxEnd = pvCol->end();
81  for (reco::VertexCollection::const_iterator vtxIter = pvCol->begin(); vtxEnd != vtxIter; ++vtxIter) {
82  if (!vtxIter->isFake() && vtxIter->ndof()>=fVtxNdofCut && std::abs(vtxIter->z())<=fVtxZCut)
83  npv++;
84  }
85 
86  //Fill the reco objects
87  fRecoObjCollection.clear();
88  fRecoObjCollection.reserve(pfCol->size());
89  for(auto const& aPF : *pfCol) {
90  RecoObj pReco;
91  pReco.pt = aPF.pt();
92  pReco.eta = aPF.eta();
93  pReco.phi = aPF.phi();
94  pReco.m = aPF.mass();
95  pReco.rapidity = aPF.rapidity();
96  pReco.charge = aPF.charge();
97  const reco::Vertex *closestVtx = nullptr;
98  double pDZ = -9999;
99  double pD0 = -9999;
100  int pVtxId = -9999;
101  bool lFirst = true;
102  const pat::PackedCandidate *lPack = dynamic_cast<const pat::PackedCandidate*>(&aPF);
103  if(lPack == nullptr ) {
104 
105  const reco::PFCandidate *pPF = dynamic_cast<const reco::PFCandidate*>(&aPF);
106  double curdz = 9999;
107  int closestVtxForUnassociateds = -9999;
108  const reco::TrackRef aTrackRef = pPF->trackRef();
109  for(auto const& aV : *pvCol) {
110  if(lFirst) {
111  if ( aTrackRef.isNonnull()) {
112  pDZ = aTrackRef->dz(aV.position());
113  pD0 = aTrackRef->d0();
114  } else if (pPF->gsfTrackRef().isNonnull()) {
115  pDZ = pPF->gsfTrackRef()->dz(aV.position());
116  pD0 = pPF->gsfTrackRef()->d0();
117  }
118  lFirst = false;
119  if(pDZ > -9999) pVtxId = 0;
120  }
121  if(aTrackRef.isNonnull() && aV.trackWeight(pPF->trackRef())>0) {
122  closestVtx = &aV;
123  break;
124  }
125  // in case it's unassocciated, keep more info
126  double tmpdz = 99999;
127  if ( aTrackRef.isNonnull() ) tmpdz = aTrackRef ->dz(aV.position());
128  else if ( pPF->gsfTrackRef().isNonnull() ) tmpdz = pPF->gsfTrackRef()->dz(aV.position());
129  if (std::abs(tmpdz) < curdz){
130  curdz = std::abs(tmpdz);
131  closestVtxForUnassociateds = pVtxId;
132  }
133  pVtxId++;
134 
135  }
136  int tmpFromPV = 0;
137  // mocking the miniAOD definitions
138  if (std::abs(pReco.charge) > 0){
139  if (closestVtx != nullptr && pVtxId > 0) tmpFromPV = 0;
140  if (closestVtx != nullptr && pVtxId == 0) tmpFromPV = 3;
141  if (closestVtx == nullptr && closestVtxForUnassociateds == 0) tmpFromPV = 2;
142  if (closestVtx == nullptr && closestVtxForUnassociateds != 0) tmpFromPV = 1;
143  }
144  pReco.dZ = pDZ;
145  pReco.d0 = pD0;
146  pReco.id = 0;
147  if (std::abs(pReco.charge) == 0){ pReco.id = 0; }
148  else{
149  if (tmpFromPV == 0){ pReco.id = 2; } // 0 is associated to PU vertex
150  else if (tmpFromPV == 3){ pReco.id = 1; }
151  else if (tmpFromPV == 1 || tmpFromPV == 2){
152  pReco.id = 0;
153  if ((fPtMaxCharged > 0) and (pReco.pt > fPtMaxCharged))
154  pReco.id = 1;
155  else if (std::abs(pReco.eta) > fEtaMaxCharged)
156  pReco.id = 1;
157  else if (fUseDZ)
158  pReco.id = (std::abs(pDZ) < fDZCut) ? 1 : 2;
159  else if (fUseFromPVLooseTight && tmpFromPV == 1)
160  pReco.id = 2;
161  else if (fUseFromPVLooseTight && tmpFromPV == 2)
162  pReco.id = 1;
163  }
164  }
165  }
166  else if(lPack->vertexRef().isNonnull() ) {
167  pDZ = lPack->dz();
168  pD0 = lPack->dxy();
169  pReco.dZ = pDZ;
170  pReco.d0 = pD0;
171 
172  pReco.id = 0;
173  if (std::abs(pReco.charge) == 0){ pReco.id = 0; }
174  if (std::abs(pReco.charge) > 0){
175  if (lPack->fromPV() == 0){ pReco.id = 2; } // 0 is associated to PU vertex
176  else if (lPack->fromPV() == (pat::PackedCandidate::PVUsedInFit)){ pReco.id = 1; }
177  else if (lPack->fromPV() == (pat::PackedCandidate::PVTight) || lPack->fromPV() == (pat::PackedCandidate::PVLoose)){
178  pReco.id = 0;
179  if ((fPtMaxCharged > 0) and (pReco.pt > fPtMaxCharged))
180  pReco.id = 1;
181  else if (std::abs(pReco.eta) > fEtaMaxCharged)
182  pReco.id = 1;
183  else if (fUseDZ)
184  pReco.id = (std::abs(pDZ) < fDZCut) ? 1 : 2;
186  pReco.id = 2;
188  pReco.id = 1;
189  }
190  }
191  }
192 
193  fRecoObjCollection.push_back(pReco);
194 
195  }
196 
197  fPuppiContainer->initialize(fRecoObjCollection);
198  fPuppiContainer->setNPV( npv );
199 
200  std::vector<double> lWeights;
201  std::vector<PuppiCandidate> lCandidates;
202  if (!fUseExistingWeights){
203  //Compute the weights and get the particles
204  lWeights = fPuppiContainer->puppiWeights();
205  lCandidates = fPuppiContainer->puppiParticles();
206  }
207  else{
208  //Use the existing weights
209  int lPackCtr = 0;
210  for(auto const& aPF : *pfCol) {
211  const pat::PackedCandidate *lPack = dynamic_cast<const pat::PackedCandidate*>(&aPF);
212  float curpupweight = -1.;
213  if(lPack == nullptr ) {
214  // throw error
215  throw edm::Exception(edm::errors::LogicError,"PuppiProducer: cannot get weights since inputs are not PackedCandidates");
216  }
217  else{
218  if (fUseWeightsNoLep){ curpupweight = lPack->puppiWeightNoLep(); }
219  else{ curpupweight = lPack->puppiWeight(); }
220  }
221  lWeights.push_back(curpupweight);
222  PuppiCandidate curjet( curpupweight*lPack->px(), curpupweight*lPack->py(), curpupweight*lPack->pz(), curpupweight*lPack->energy());
223  curjet.set_user_index(lPackCtr);
224  lCandidates.push_back(curjet);
225  lPackCtr++;
226  }
227  }
228 
229  //Fill it into the event
230  std::unique_ptr<edm::ValueMap<float> > lPupOut(new edm::ValueMap<float>());
231  edm::ValueMap<float>::Filler lPupFiller(*lPupOut);
232  lPupFiller.insert(hPFProduct,lWeights.begin(),lWeights.end());
233  lPupFiller.fill();
234 
235  // This is a dummy to access the "translate" method which is a
236  // non-static member function even though it doesn't need to be.
237  // Will fix in the future.
238  static const reco::PFCandidate dummySinceTranslateIsNotStatic;
239 
240  // Fill a new PF/Packed Candidate Collection and write out the ValueMap of the new p4s.
241  // Since the size of the ValueMap must be equal to the input collection, we need
242  // to search the "puppi" particles to find a match for each input. If none is found,
243  // the input is set to have a four-vector of 0,0,0,0
246  std::unique_ptr<edm::ValueMap<LorentzVector> > p4PupOut(new edm::ValueMap<LorentzVector>());
247  LorentzVectorCollection puppiP4s;
248  std::vector<reco::CandidatePtr> values(hPFProduct->size());
249 
250  int val = -1;
251  for ( auto const& aCand : *hPFProduct) {
252  val++;
253  std::unique_ptr<pat::PackedCandidate> pCand;
254  std::unique_ptr<reco::PFCandidate> pfCand;
256  const pat::PackedCandidate *cand = dynamic_cast<const pat::PackedCandidate*>(&aCand);
257  if(!cand)
258  throw edm::Exception(edm::errors::LogicError,"PuppiProducer: inputs are not PackedCandidates");
259  pCand.reset( new pat::PackedCandidate(*cand) );
260  } else {
261  auto id = dummySinceTranslateIsNotStatic.translatePdgIdToType(aCand.pdgId());
262  const reco::PFCandidate *cand = dynamic_cast<const reco::PFCandidate*>(&aCand);
263  pfCand.reset( new reco::PFCandidate( cand ? *cand : reco::PFCandidate(aCand.charge(), aCand.p4(), id) ) );
264  }
265  LorentzVector pVec;
266 
267  //get an index to a pup in lCandidates: either fUseExistingWeights with no skips or get from fPuppiContainer
268  int iPuppiMatched = fUseExistingWeights ? val : fPuppiContainer->recoToPup()[val];
269  if ( iPuppiMatched >= 0 ) {
270  auto const& puppiMatched = lCandidates[iPuppiMatched];
271  pVec.SetPxPyPzE(puppiMatched.px(),puppiMatched.py(),puppiMatched.pz(),puppiMatched.E());
273  if(fPuppiForLeptons)
274  pCand->setPuppiWeight(pCand->puppiWeight(),lWeights[val]);
275  else
276  pCand->setPuppiWeight(lWeights[val],pCand->puppiWeightNoLep());
277  }
278  } else {
279  pVec.SetPxPyPzE( 0, 0, 0, 0);
281  pCand->setPuppiWeight(0,0);
282  }
283  }
284  puppiP4s.push_back( pVec );
285 
287  pCand->setP4(pVec);
288  pCand->setSourceCandidatePtr( aCand.sourceCandidatePtr(0) );
289  fPackedPuppiCandidates->push_back(*pCand);
290  } else {
291  pfCand->setP4(pVec);
292  pfCand->setSourceCandidatePtr( aCand.sourceCandidatePtr(0) );
293  fPuppiCandidates->push_back(*pfCand);
294  }
295  }
296 
297  //Compute the modified p4s
298  edm::ValueMap<LorentzVector>::Filler p4PupFiller(*p4PupOut);
299  p4PupFiller.insert(hPFProduct,puppiP4s.begin(), puppiP4s.end() );
300  p4PupFiller.fill();
301 
302  iEvent.put(std::move(lPupOut));
303  iEvent.put(std::move(p4PupOut));
306  for(unsigned int ic=0, nc = oh->size(); ic < nc; ++ic) {
307  reco::CandidatePtr pkref( oh, ic );
308  values[ic] = pkref;
309  }
310  } else {
312  for(unsigned int ic=0, nc = oh->size(); ic < nc; ++ic) {
313  reco::CandidatePtr pkref( oh, ic );
314  values[ic] = pkref;
315  }
316  }
317  std::unique_ptr<edm::ValueMap<reco::CandidatePtr> > pfMap_p(new edm::ValueMap<reco::CandidatePtr>());
319  filler.insert(hPFProduct, values.begin(), values.end());
320  filler.fill();
321  iEvent.put(std::move(pfMap_p));
322 
323 
326 
327  // all the different alphas per particle
328  // THE alpha per particle
329  std::unique_ptr<std::vector<double> > theAlphas(new std::vector<double>(fPuppiContainer->puppiAlphas()));
330  std::unique_ptr<std::vector<double> > theAlphasMed(new std::vector<double>(fPuppiContainer->puppiAlphasMed()));
331  std::unique_ptr<std::vector<double> > theAlphasRms(new std::vector<double>(fPuppiContainer->puppiAlphasRMS()));
332  std::unique_ptr<std::vector<double> > alphas(new std::vector<double>(fPuppiContainer->puppiRawAlphas()));
333  std::unique_ptr<double> nalgos(new double(fPuppiContainer->puppiNAlgos()));
334 
335  iEvent.put(std::move(alphas),"PuppiRawAlphas");
336  iEvent.put(std::move(nalgos),"PuppiNAlgos");
337  iEvent.put(std::move(theAlphas),"PuppiAlphas");
338  iEvent.put(std::move(theAlphasMed),"PuppiAlphasMed");
339  iEvent.put(std::move(theAlphasRms),"PuppiAlphasRms");
340  }
341 
342 }
343 
344 // ------------------------------------------------------------------------------------------
346 }
347 // ------------------------------------------------------------------------------------------
349 }
350 // ------------------------------------------------------------------------------------------
353  desc.add<bool>("puppiDiagnostics", false);
354  desc.add<bool>("puppiForLeptons", false);
355  desc.add<bool>("UseFromPVLooseTight", false);
356  desc.add<bool>("UseDeltaZCut", true);
357  desc.add<double>("DeltaZCut", 0.3);
358  desc.add<double>("PtMaxCharged", 0.);
359  desc.add<double>("EtaMaxCharged", 99999.);
360  desc.add<double>("PtMaxNeutrals", 200.);
361  desc.add<double>("PtMaxNeutralsStartSlope", 0.);
362  desc.add<bool>("useExistingWeights", false);
363  desc.add<bool>("useWeightsNoLep", false);
364  desc.add<bool>("clonePackedCands", false);
365  desc.add<int>("vtxNdofCut", 4);
366  desc.add<double>("vtxZCut", 24);
367  desc.add<edm::InputTag>("candName", edm::InputTag("particleFlow"));
368  desc.add<edm::InputTag>("vertexName", edm::InputTag("offlinePrimaryVertices"));
369  desc.add<bool>("applyCHS", true);
370  desc.add<bool>("invertPuppi", false);
371  desc.add<bool>("useExp", false);
372  desc.add<double>("MinPuppiWeight", .01);
373 
375 
376  descriptions.add("PuppiProducer", desc);
377 }
378 //define this as a plug-in
float puppiWeight() const
virtual float dz(size_t ipv=0) const
dz with respect to the PV[ipv]
T getParameter(std::string const &) const
bool fUseFromPVLooseTight
Definition: PuppiProducer.h:51
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:251
double px() const override
x coordinate of momentum vector
float puppiWeightNoLep() const
Weight from full PUPPI.
float d0
Definition: RecoObj.h:28
virtual void endJob()
bool fUseWeightsNoLep
Definition: PuppiProducer.h:57
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
std::vector< reco::PFCandidate > PFOutputCollection
Definition: PuppiProducer.h:35
static void fillDescriptionsPuppiAlgo(edm::ParameterSetDescription &desc)
Definition: PuppiAlgo.cc:217
Definition: RecoObj.h:4
bool fPuppiDiagnostics
Definition: PuppiProducer.h:49
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:53
std::vector< pat::PackedCandidate > PackedOutputCollection
Definition: PuppiProducer.h:36
std::vector< RecoObj > fRecoObjCollection
Definition: PuppiProducer.h:62
size_type size() const
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
std::vector< LorentzVector > LorentzVectorCollection
Definition: PuppiProducer.h:31
virtual void beginJob()
math::XYZTLorentzVector LorentzVector
Definition: PuppiProducer.h:30
const reco::VertexRef vertexRef() const
double fEtaMaxCharged
Definition: PuppiProducer.h:55
std::unique_ptr< PuppiContainer > fPuppiContainer
Definition: PuppiProducer.h:61
int charge
Definition: RecoObj.h:29
double pz() const override
z coordinate of momentum vector
edm::EDGetTokenT< VertexCollection > tokenVertices_
Definition: PuppiProducer.h:45
float dZ
Definition: RecoObj.h:27
double py() const override
y coordinate of momentum vector
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:442
float rapidity
Definition: RecoObj.h:17
float pt
Definition: RecoObj.h:17
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
int id
Definition: RecoObj.h:18
float phi
Definition: RecoObj.h:17
const PVAssoc fromPV(size_t ipv=0) const
std::unique_ptr< PackedOutputCollection > fPackedPuppiCandidates
Definition: PuppiProducer.h:64
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool fClonePackedCands
Definition: PuppiProducer.h:58
double energy() const override
energy
void produce(edm::Event &, const edm::EventSetup &) override
bool fPuppiForLeptons
Definition: PuppiProducer.h:50
ParameterDescriptionBase * add(U const &iLabel, T const &value)
PuppiProducer(const edm::ParameterSet &)
float m
Definition: RecoObj.h:17
T const * product() const
Definition: Handle.h:74
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
double fPtMaxCharged
Definition: PuppiProducer.h:54
void add(std::string const &label, ParameterSetDescription const &psetDescription)
ParticleType translatePdgIdToType(int pdgid) const
Definition: PFCandidate.cc:233
float eta
Definition: RecoObj.h:17
bool fUseExistingWeights
Definition: PuppiProducer.h:56
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:40
reco::GsfTrackRef gsfTrackRef() const
Definition: PFCandidate.cc:480
edm::EDGetTokenT< CandidateView > tokenPFCandidates_
Definition: PuppiProducer.h:44
virtual float dxy() const
dxy with respect to the PV ref
~PuppiProducer() override
std::unique_ptr< PFOutputCollection > fPuppiCandidates
Definition: PuppiProducer.h:63
def move(src, dest)
Definition: eostools.py:511