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