CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PuppiProducer.cc
Go to the documentation of this file.
1 // system include files
2 #include <memory>
3 
4 // user include files
24 //Main File
25 #include "fastjet/PseudoJet.hh"
27 
28 
29 
30 // ------------------------------------------------------------------------------------------
32  fPuppiDiagnostics = iConfig.getParameter<bool>("puppiDiagnostics");
33  fPuppiForLeptons = iConfig.getParameter<bool>("puppiForLeptons");
34  fUseDZ = iConfig.getParameter<bool>("UseDeltaZCut");
35  fDZCut = iConfig.getParameter<double>("DeltaZCut");
36  fUseExistingWeights = iConfig.getParameter<bool>("useExistingWeights");
37  fUseWeightsNoLep = iConfig.getParameter<bool>("useWeightsNoLep");
38  fClonePackedCands = iConfig.getParameter<bool>("clonePackedCands");
39  fVtxNdofCut = iConfig.getParameter<int>("vtxNdofCut");
40  fVtxZCut = iConfig.getParameter<double>("vtxZCut");
41  fPuppiContainer = std::unique_ptr<PuppiContainer> ( new PuppiContainer(iConfig) );
42 
44  = consumes<CandidateView>(iConfig.getParameter<edm::InputTag>("candName"));
46  = consumes<VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexName"));
47 
48 
49  produces<edm::ValueMap<float> > ();
50  produces<edm::ValueMap<LorentzVector> > ();
51  produces< edm::ValueMap<reco::CandidatePtr> >();
52 
53  if (fUseExistingWeights || fClonePackedCands)
54  produces<pat::PackedCandidateCollection>();
55  else
56  produces<reco::PFCandidateCollection>();
57 
58  if (fPuppiDiagnostics){
59  produces<double> ("PuppiNAlgos");
60  produces<std::vector<double>> ("PuppiRawAlphas");
61  produces<std::vector<double>> ("PuppiAlphas");
62  produces<std::vector<double>> ("PuppiAlphasMed");
63  produces<std::vector<double>> ("PuppiAlphasRms");
64  }
65 }
66 // ------------------------------------------------------------------------------------------
68 }
69 // ------------------------------------------------------------------------------------------
71 
72  // Get PFCandidate Collection
73  edm::Handle<CandidateView> hPFProduct;
74  iEvent.getByToken(tokenPFCandidates_,hPFProduct);
75  const CandidateView *pfCol = hPFProduct.product();
76 
77  // Get vertex collection w/PV as the first entry?
79  iEvent.getByToken(tokenVertices_,hVertexProduct);
80  const reco::VertexCollection *pvCol = hVertexProduct.product();
81 
82  int npv = 0;
83  const reco::VertexCollection::const_iterator vtxEnd = pvCol->end();
84  for (reco::VertexCollection::const_iterator vtxIter = pvCol->begin(); vtxEnd != vtxIter; ++vtxIter) {
85  if (!vtxIter->isFake() && vtxIter->ndof()>=fVtxNdofCut && fabs(vtxIter->z())<=fVtxZCut)
86  npv++;
87  }
88 
89  //Fill the reco objects
90  fRecoObjCollection.clear();
91  for(CandidateView::const_iterator itPF = pfCol->begin(); itPF!=pfCol->end(); itPF++) {
92  // std::cout << "itPF->pdgId() = " << itPF->pdgId() << std::endl;
93  RecoObj pReco;
94  pReco.pt = itPF->pt();
95  pReco.eta = itPF->eta();
96  pReco.phi = itPF->phi();
97  pReco.m = itPF->mass();
98  pReco.rapidity = itPF->rapidity();
99  pReco.charge = itPF->charge();
100  const reco::Vertex *closestVtx = 0;
101  double pDZ = -9999;
102  double pD0 = -9999;
103  int pVtxId = -9999;
104  bool lFirst = true;
105  const pat::PackedCandidate *lPack = dynamic_cast<const pat::PackedCandidate*>(&(*itPF));
106  if(lPack == 0 ) {
107 
108  const reco::PFCandidate *pPF = dynamic_cast<const reco::PFCandidate*>(&(*itPF));
109  double curdz = 9999;
110  int closestVtxForUnassociateds = -9999;
111  for(reco::VertexCollection::const_iterator iV = pvCol->begin(); iV!=pvCol->end(); ++iV) {
112  if(lFirst) {
113  if ( pPF->trackRef().isNonnull() ) pDZ = pPF->trackRef() ->dz(iV->position());
114  else if ( pPF->gsfTrackRef().isNonnull() ) pDZ = pPF->gsfTrackRef()->dz(iV->position());
115  if ( pPF->trackRef().isNonnull() ) pD0 = pPF->trackRef() ->d0();
116  else if ( pPF->gsfTrackRef().isNonnull() ) pD0 = pPF->gsfTrackRef()->d0();
117  lFirst = false;
118  if(pDZ > -9999) pVtxId = 0;
119  }
120  if(iV->trackWeight(pPF->trackRef())>0) {
121  closestVtx = &(*iV);
122  break;
123  }
124  // in case it's unassocciated, keep more info
125  double tmpdz = 99999;
126  if ( pPF->trackRef().isNonnull() ) tmpdz = pPF->trackRef() ->dz(iV->position());
127  else if ( pPF->gsfTrackRef().isNonnull() ) tmpdz = pPF->gsfTrackRef()->dz(iV->position());
128  if (fabs(tmpdz) < curdz){
129  curdz = fabs(tmpdz);
130  closestVtxForUnassociateds = pVtxId;
131  }
132  pVtxId++;
133 
134  }
135  int tmpFromPV = 0;
136  // mocking the miniAOD definitions
137  if (closestVtx != 0 && fabs(pReco.charge) > 0 && pVtxId > 0) tmpFromPV = 0;
138  if (closestVtx != 0 && fabs(pReco.charge) > 0 && pVtxId == 0) tmpFromPV = 3;
139  if (closestVtx == 0 && fabs(pReco.charge) > 0 && closestVtxForUnassociateds == 0) tmpFromPV = 2;
140  if (closestVtx == 0 && fabs(pReco.charge) > 0 && closestVtxForUnassociateds != 0) tmpFromPV = 1;
141  pReco.dZ = pDZ;
142  pReco.d0 = pD0;
143  pReco.id = 0;
144  if (fabs(pReco.charge) == 0){ pReco.id = 0; }
145  else{
146  if (tmpFromPV == 0){ pReco.id = 2; } // 0 is associated to PU vertex
147  if (tmpFromPV == 3){ pReco.id = 1; }
148  if (tmpFromPV == 1 || tmpFromPV == 2){
149  pReco.id = 0;
150  if (!fPuppiForLeptons && fUseDZ && (fabs(pDZ) < fDZCut)) pReco.id = 1;
151  if (!fPuppiForLeptons && fUseDZ && (fabs(pDZ) > fDZCut)) pReco.id = 2;
152  if (fPuppiForLeptons && tmpFromPV == 1) pReco.id = 2;
153  if (fPuppiForLeptons && tmpFromPV == 2) pReco.id = 1;
154  }
155  }
156  }
157  else if(lPack->vertexRef().isNonnull() ) {
158  pDZ = lPack->dz();
159  pD0 = lPack->dxy();
160  closestVtx = &(*(lPack->vertexRef()));
161  pReco.dZ = pDZ;
162  pReco.d0 = pD0;
163 
164  pReco.id = 0;
165  if (fabs(pReco.charge) == 0){ pReco.id = 0; }
166  if (fabs(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 && (fabs(pDZ) < fDZCut)) pReco.id = 1;
172  if (!fPuppiForLeptons && fUseDZ && (fabs(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<fastjet::PseudoJet> 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(CandidateView::const_iterator itPF = pfCol->begin(); itPF!=pfCol->end(); itPF++) {
197  const pat::PackedCandidate *lPack = dynamic_cast<const pat::PackedCandidate*>(&(*itPF));
198  float curpupweight = -1.;
199  if(lPack == 0 ) {
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 = itPF->puppiWeightNoLep(); }
205  // else{ curpupweight = itPF->puppiWeight(); }
206  curpupweight = lPack->puppiWeight();
207  }
208  lWeights.push_back(curpupweight);
209  fastjet::PseudoJet curjet( curpupweight*lPack->px(), curpupweight*lPack->py(), curpupweight*lPack->pz(), curpupweight*lPack->energy());
210  curjet.set_user_index(lPackCtr);
211  lCandidates.push_back(curjet);
212  lPackCtr++;
213  }
214  }
215 
216  //Fill it into the event
217  std::auto_ptr<edm::ValueMap<float> > lPupOut(new edm::ValueMap<float>());
218  edm::ValueMap<float>::Filler lPupFiller(*lPupOut);
219  lPupFiller.insert(hPFProduct,lWeights.begin(),lWeights.end());
220  lPupFiller.fill();
221 
222  // This is a dummy to access the "translate" method which is a
223  // non-static member function even though it doesn't need to be.
224  // Will fix in the future.
225  static const reco::PFCandidate dummySinceTranslateIsNotStatic;
226 
227  // Fill a new PF/Packed Candidate Collection and write out the ValueMap of the new p4s.
228  // Since the size of the ValueMap must be equal to the input collection, we need
229  // to search the "puppi" particles to find a match for each input. If none is found,
230  // the input is set to have a four-vector of 0,0,0,0
233  std::auto_ptr<edm::ValueMap<LorentzVector> > p4PupOut(new edm::ValueMap<LorentzVector>());
234  LorentzVectorCollection puppiP4s;
235  std::vector<reco::CandidatePtr> values(hPFProduct->size());
236 
237  for ( auto i0 = hPFProduct->begin(),
238  i0begin = hPFProduct->begin(),
239  i0end = hPFProduct->end(); i0 != i0end; ++i0 ) {
240  std::unique_ptr<pat::PackedCandidate> pCand;
241  std::unique_ptr<reco::PFCandidate> pfCand;
243  const pat::PackedCandidate *cand = dynamic_cast<const pat::PackedCandidate*>(&(*i0));
244  if(!cand)
245  throw edm::Exception(edm::errors::LogicError,"PuppiProducer: inputs are not PackedCandidates");
246  pCand.reset( new pat::PackedCandidate(*cand) );
247  } else {
248  auto id = dummySinceTranslateIsNotStatic.translatePdgIdToType(i0->pdgId());
249  const reco::PFCandidate *cand = dynamic_cast<const reco::PFCandidate*>(&(*i0));
250  pfCand.reset( new reco::PFCandidate( cand ? *cand : reco::PFCandidate(i0->charge(), i0->p4(), id) ) );
251  }
252  LorentzVector pVec = i0->p4();
253  int val = i0 - i0begin;
254 
255  // Find the Puppi particle matched to the input collection using the "user_index" of the object.
256  auto puppiMatched = find_if( lCandidates.begin(), lCandidates.end(), [&val]( fastjet::PseudoJet const & i ){ return i.user_index() == val; } );
257  if ( puppiMatched != lCandidates.end() ) {
258  pVec.SetPxPyPzE(puppiMatched->px(),puppiMatched->py(),puppiMatched->pz(),puppiMatched->E());
259  } else {
260  pVec.SetPxPyPzE( 0, 0, 0, 0);
261  }
262  puppiP4s.push_back( pVec );
263 
265  pCand->setP4(pVec);
266  pCand->setSourceCandidatePtr( i0->sourceCandidatePtr(0) );
267  fPackedPuppiCandidates->push_back(*pCand);
268  } else {
269  pfCand->setP4(pVec);
270  pfCand->setSourceCandidatePtr( i0->sourceCandidatePtr(0) );
271  fPuppiCandidates->push_back(*pfCand);
272  }
273  }
274 
275  //Compute the modified p4s
276  edm::ValueMap<LorentzVector>::Filler p4PupFiller(*p4PupOut);
277  p4PupFiller.insert(hPFProduct,puppiP4s.begin(), puppiP4s.end() );
278  p4PupFiller.fill();
279 
280  iEvent.put(lPupOut);
281  iEvent.put(p4PupOut);
284  for(unsigned int ic=0, nc = oh->size(); ic < nc; ++ic) {
285  reco::CandidatePtr pkref( oh, ic );
286  values[ic] = pkref;
287  }
288  } else {
290  for(unsigned int ic=0, nc = oh->size(); ic < nc; ++ic) {
291  reco::CandidatePtr pkref( oh, ic );
292  values[ic] = pkref;
293  }
294  }
295  std::auto_ptr<edm::ValueMap<reco::CandidatePtr> > pfMap_p(new edm::ValueMap<reco::CandidatePtr>());
297  filler.insert(hPFProduct, values.begin(), values.end());
298  filler.fill();
299  iEvent.put(pfMap_p);
300 
301 
304 
305  // all the different alphas per particle
306  // THE alpha per particle
307  std::auto_ptr<std::vector<double> > theAlphas(new std::vector<double>(fPuppiContainer->puppiAlphas()));
308  std::auto_ptr<std::vector<double> > theAlphasMed(new std::vector<double>(fPuppiContainer->puppiAlphasMed()));
309  std::auto_ptr<std::vector<double> > theAlphasRms(new std::vector<double>(fPuppiContainer->puppiAlphasRMS()));
310  std::auto_ptr<std::vector<double> > alphas(new std::vector<double>(fPuppiContainer->puppiRawAlphas()));
311  std::auto_ptr<double> nalgos(new double(fPuppiContainer->puppiNAlgos()));
312 
313  iEvent.put(alphas,"PuppiRawAlphas");
314  iEvent.put(nalgos,"PuppiNAlgos");
315  iEvent.put(theAlphas,"PuppiAlphas");
316  iEvent.put(theAlphasMed,"PuppiAlphasMed");
317  iEvent.put(theAlphasRms,"PuppiAlphasRms");
318  }
319 
320 }
321 
322 // ------------------------------------------------------------------------------------------
324 }
325 // ------------------------------------------------------------------------------------------
327 }
328 // ------------------------------------------------------------------------------------------
330 }
331 // ------------------------------------------------------------------------------------------
333 }
334 // ------------------------------------------------------------------------------------------
336 }
337 // ------------------------------------------------------------------------------------------
339 }
340 // ------------------------------------------------------------------------------------------
342  //The following says we do not know what parameters are allowed so do no validation
343  // Please change this to state exactly what you do use, even if it is no parameters
345  desc.setUnknown();
346  descriptions.addDefault(desc);
347 }
348 //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
int i
Definition: DBlmapReader.cc:9
virtual void produce(edm::Event &, const edm::EventSetup &)
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
float d0
Definition: RecoObj.h:28
virtual void endJob()
std::auto_ptr< PackedOutputCollection > fPackedPuppiCandidates
Definition: PuppiProducer.h:62
bool fUseWeightsNoLep
Definition: PuppiProducer.h:55
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
std::vector< reco::PFCandidate > PFOutputCollection
Definition: PuppiProducer.h:31
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Definition: RecoObj.h:4
bool fPuppiDiagnostics
Definition: PuppiProducer.h:50
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:52
std::vector< pat::PackedCandidate > PackedOutputCollection
Definition: PuppiProducer.h:32
std::vector< RecoObj > fRecoObjCollection
Definition: PuppiProducer.h:60
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
std::vector< LorentzVector > LorentzVectorCollection
Definition: PuppiProducer.h:27
virtual void beginJob()
math::XYZTLorentzVector LorentzVector
Definition: PuppiProducer.h:26
const reco::VertexRef vertexRef() const
std::unique_ptr< PuppiContainer > fPuppiContainer
Definition: PuppiProducer.h:59
int charge
Definition: RecoObj.h:29
edm::EDGetTokenT< VertexCollection > tokenVertices_
Definition: PuppiProducer.h:46
float dZ
Definition: RecoObj.h:27
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:433
float rapidity
Definition: RecoObj.h:17
float pt
Definition: RecoObj.h:17
int iEvent
Definition: GenABIO.cc:230
void addDefault(ParameterSetDescription const &psetDescription)
const_iterator begin() const
int id
Definition: RecoObj.h:18
float phi
Definition: RecoObj.h:17
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
const PVAssoc fromPV(size_t ipv=0) const
virtual double py() const
y coordinate of momentum vector
virtual double energy() const
energy
bool fClonePackedCands
Definition: PuppiProducer.h:56
bool fPuppiForLeptons
Definition: PuppiProducer.h:51
PuppiProducer(const edm::ParameterSet &)
virtual void beginRun(edm::Run &, edm::EventSetup const &)
virtual void endLuminosityBlock(edm::LuminosityBlock &, edm::EventSetup const &)
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:224
float eta
Definition: RecoObj.h:17
bool fUseExistingWeights
Definition: PuppiProducer.h:54
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
reco::GsfTrackRef gsfTrackRef() const
Definition: PFCandidate.cc:471
edm::EDGetTokenT< CandidateView > tokenPFCandidates_
Definition: PuppiProducer.h:45
std::auto_ptr< PFOutputCollection > fPuppiCandidates
Definition: PuppiProducer.h:61
const_iterator end() const
virtual void endRun(edm::Run &, edm::EventSetup const &)
virtual double px() const
x coordinate of momentum vector
virtual void beginLuminosityBlock(edm::LuminosityBlock &, edm::EventSetup const &)
virtual float dxy() const
dxy with respect to the PV ref
virtual double pz() const
z coordinate of momentum vector
virtual const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
Definition: Run.h:43