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  fVtxNdofCut = iConfig.getParameter<int>("vtxNdofCut");
39  fVtxZCut = iConfig.getParameter<double>("vtxZCut");
40  fPuppiContainer = std::unique_ptr<PuppiContainer> ( new PuppiContainer(iConfig) );
41 
43  = consumes<CandidateView>(iConfig.getParameter<edm::InputTag>("candName"));
45  = consumes<VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexName"));
46 
47 
48  produces<edm::ValueMap<float> > ();
49  produces<edm::ValueMap<LorentzVector> > ();
50  produces< edm::ValueMap<reco::CandidatePtr> >();
51 
52  produces<PFOutputCollection>();
53 
54  if (fPuppiDiagnostics){
55  produces<double> ("PuppiNAlgos");
56  produces<std::vector<double>> ("PuppiRawAlphas");
57  produces<std::vector<double>> ("PuppiAlphas");
58  produces<std::vector<double>> ("PuppiAlphasMed");
59  produces<std::vector<double>> ("PuppiAlphasRms");
60  }
61 }
62 // ------------------------------------------------------------------------------------------
64 }
65 // ------------------------------------------------------------------------------------------
67 
68  // Get PFCandidate Collection
69  edm::Handle<CandidateView> hPFProduct;
70  iEvent.getByToken(tokenPFCandidates_,hPFProduct);
71  const CandidateView *pfCol = hPFProduct.product();
72 
73  // Get vertex collection w/PV as the first entry?
75  iEvent.getByToken(tokenVertices_,hVertexProduct);
76  const reco::VertexCollection *pvCol = hVertexProduct.product();
77 
78  int npv = 0;
79  const reco::VertexCollection::const_iterator vtxEnd = pvCol->end();
80  for (reco::VertexCollection::const_iterator vtxIter = pvCol->begin(); vtxEnd != vtxIter; ++vtxIter) {
81  if (!vtxIter->isFake() && vtxIter->ndof()>=fVtxNdofCut && fabs(vtxIter->z())<=fVtxZCut)
82  npv++;
83  }
84 
85  //Fill the reco objects
86  fRecoObjCollection.clear();
87  for(CandidateView::const_iterator itPF = pfCol->begin(); itPF!=pfCol->end(); itPF++) {
88  // std::cout << "itPF->pdgId() = " << itPF->pdgId() << std::endl;
89  RecoObj pReco;
90  pReco.pt = itPF->pt();
91  pReco.eta = itPF->eta();
92  pReco.phi = itPF->phi();
93  pReco.m = itPF->mass();
94  pReco.rapidity = itPF->rapidity();
95  pReco.charge = itPF->charge();
96  const reco::Vertex *closestVtx = 0;
97  double pDZ = -9999;
98  double pD0 = -9999;
99  int pVtxId = -9999;
100  bool lFirst = true;
101  const pat::PackedCandidate *lPack = dynamic_cast<const pat::PackedCandidate*>(&(*itPF));
102  if(lPack == 0 ) {
103 
104  const reco::PFCandidate *pPF = dynamic_cast<const reco::PFCandidate*>(&(*itPF));
105  double curdz = 9999;
106  int closestVtxForUnassociateds = -9999;
107  for(reco::VertexCollection::const_iterator iV = pvCol->begin(); iV!=pvCol->end(); ++iV) {
108  if(lFirst) {
109  if ( pPF->trackRef().isNonnull() ) pDZ = pPF->trackRef() ->dz(iV->position());
110  else if ( pPF->gsfTrackRef().isNonnull() ) pDZ = pPF->gsfTrackRef()->dz(iV->position());
111  if ( pPF->trackRef().isNonnull() ) pD0 = pPF->trackRef() ->d0();
112  else if ( pPF->gsfTrackRef().isNonnull() ) pD0 = pPF->gsfTrackRef()->d0();
113  lFirst = false;
114  if(pDZ > -9999) pVtxId = 0;
115  }
116  if(iV->trackWeight(pPF->trackRef())>0) {
117  closestVtx = &(*iV);
118  break;
119  }
120  // in case it's unassocciated, keep more info
121  double tmpdz = 99999;
122  if ( pPF->trackRef().isNonnull() ) tmpdz = pPF->trackRef() ->dz(iV->position());
123  else if ( pPF->gsfTrackRef().isNonnull() ) tmpdz = pPF->gsfTrackRef()->dz(iV->position());
124  if (fabs(tmpdz) < curdz){
125  curdz = fabs(tmpdz);
126  closestVtxForUnassociateds = pVtxId;
127  }
128  pVtxId++;
129 
130  }
131  int tmpFromPV = 0;
132  // mocking the miniAOD definitions
133  if (closestVtx != 0 && fabs(pReco.charge) > 0 && pVtxId > 0) tmpFromPV = 0;
134  if (closestVtx != 0 && fabs(pReco.charge) > 0 && pVtxId == 0) tmpFromPV = 3;
135  if (closestVtx == 0 && fabs(pReco.charge) > 0 && closestVtxForUnassociateds == 0) tmpFromPV = 2;
136  if (closestVtx == 0 && fabs(pReco.charge) > 0 && closestVtxForUnassociateds != 0) tmpFromPV = 1;
137  pReco.dZ = pDZ;
138  pReco.d0 = pD0;
139  pReco.id = 0;
140  if (fabs(pReco.charge) == 0){ pReco.id = 0; }
141  else{
142  if (tmpFromPV == 0){ pReco.id = 2; } // 0 is associated to PU vertex
143  if (tmpFromPV == 3){ pReco.id = 1; }
144  if (tmpFromPV == 1 || tmpFromPV == 2){
145  pReco.id = 0;
146  if (!fPuppiForLeptons && fUseDZ && (fabs(pDZ) < fDZCut)) pReco.id = 1;
147  if (!fPuppiForLeptons && fUseDZ && (fabs(pDZ) > fDZCut)) pReco.id = 2;
148  if (fPuppiForLeptons && tmpFromPV == 1) pReco.id = 2;
149  if (fPuppiForLeptons && tmpFromPV == 2) pReco.id = 1;
150  }
151  }
152  }
153  else if(lPack->vertexRef().isNonnull() ) {
154  pDZ = lPack->dz();
155  pD0 = lPack->dxy();
156  closestVtx = &(*(lPack->vertexRef()));
157  pReco.dZ = pDZ;
158  pReco.d0 = pD0;
159 
160  pReco.id = 0;
161  if (fabs(pReco.charge) == 0){ pReco.id = 0; }
162  if (fabs(pReco.charge) > 0){
163  if (lPack->fromPV() == 0){ pReco.id = 2; } // 0 is associated to PU vertex
164  if (lPack->fromPV() == (pat::PackedCandidate::PVUsedInFit)){ pReco.id = 1; }
165  if (lPack->fromPV() == (pat::PackedCandidate::PVTight) || lPack->fromPV() == (pat::PackedCandidate::PVLoose)){
166  pReco.id = 0;
167  if (!fPuppiForLeptons && fUseDZ && (fabs(pDZ) < fDZCut)) pReco.id = 1;
168  if (!fPuppiForLeptons && fUseDZ && (fabs(pDZ) > fDZCut)) pReco.id = 2;
169  if (fPuppiForLeptons && lPack->fromPV() == (pat::PackedCandidate::PVLoose)) pReco.id = 2;
170  if (fPuppiForLeptons && lPack->fromPV() == (pat::PackedCandidate::PVTight)) pReco.id = 1;
171  }
172  }
173  }
174 
175  fRecoObjCollection.push_back(pReco);
176 
177  }
178 
179  fPuppiContainer->initialize(fRecoObjCollection);
180  fPuppiContainer->setNPV( npv );
181 
182  std::vector<double> lWeights;
183  std::vector<fastjet::PseudoJet> lCandidates;
184  if (!fUseExistingWeights){
185  //Compute the weights and get the particles
186  lWeights = fPuppiContainer->puppiWeights();
187  lCandidates = fPuppiContainer->puppiParticles();
188  }
189  else{
190  //Use the existing weights
191  int lPackCtr = 0;
192  for(CandidateView::const_iterator itPF = pfCol->begin(); itPF!=pfCol->end(); itPF++) {
193  const pat::PackedCandidate *lPack = dynamic_cast<const pat::PackedCandidate*>(&(*itPF));
194  float curpupweight = -1.;
195  if(lPack == 0 ) {
196  // throw error
197  throw edm::Exception(edm::errors::LogicError,"PuppiProducer: cannot get weights since inputs are not packedPFCandidates");
198  }
199  else{
200  // if (fUseWeightsNoLep){ curpupweight = itPF->puppiWeightNoLep(); }
201  // else{ curpupweight = itPF->puppiWeight(); }
202  curpupweight = lPack->puppiWeight();
203  }
204  lWeights.push_back(curpupweight);
205  fastjet::PseudoJet curjet( curpupweight*lPack->px(), curpupweight*lPack->py(), curpupweight*lPack->pz(), curpupweight*lPack->energy());
206  curjet.set_user_index(lPackCtr);
207  lCandidates.push_back(curjet);
208  lPackCtr++;
209  }
210  }
211 
212  //Fill it into the event
213  std::auto_ptr<edm::ValueMap<float> > lPupOut(new edm::ValueMap<float>());
214  edm::ValueMap<float>::Filler lPupFiller(*lPupOut);
215  lPupFiller.insert(hPFProduct,lWeights.begin(),lWeights.end());
216  lPupFiller.fill();
217 
218  // This is a dummy to access the "translate" method which is a
219  // non-static member function even though it doesn't need to be.
220  // Will fix in the future.
221  static const reco::PFCandidate dummySinceTranslateIsNotStatic;
222 
223  //Fill a new PF Candidate Collection and write out the ValueMap of the new p4s.
224  // Since the size of the ValueMap must be equal to the input collection, we need
225  // to search the "puppi" particles to find a match for each input. If none is found,
226  // the input is set to have a four-vector of 0,0,0,0
228  std::auto_ptr<edm::ValueMap<LorentzVector> > p4PupOut(new edm::ValueMap<LorentzVector>());
229  LorentzVectorCollection puppiP4s;
230  std::vector<reco::CandidatePtr> values(hPFProduct->size());
231  //std::vector<int> values(hPFProduct->size());
232 
233  for ( auto i0 = hPFProduct->begin(),
234  i0begin = hPFProduct->begin(),
235  i0end = hPFProduct->end(); i0 != i0end; ++i0 ) {
236  //for(unsigned int i0 = 0; i0 < lCandidates.size(); i0++) {
237  //reco::PFCandidate pCand;
238  auto id = dummySinceTranslateIsNotStatic.translatePdgIdToType(i0->pdgId());
239  const reco::PFCandidate *pPF = dynamic_cast<const reco::PFCandidate*>(&(*i0));
240  reco::PFCandidate pCand( pPF ? *pPF : reco::PFCandidate(i0->charge(), i0->p4(), id) );
241  LorentzVector pVec = i0->p4();
242  int val = i0 - i0begin;
243 
244  // Find the Puppi particle matched to the input collection using the "user_index" of the object.
245  auto puppiMatched = find_if( lCandidates.begin(), lCandidates.end(), [&val]( fastjet::PseudoJet const & i ){ return i.user_index() == val; } );
246  if ( puppiMatched != lCandidates.end() ) {
247  pVec.SetPxPyPzE(puppiMatched->px(),puppiMatched->py(),puppiMatched->pz(),puppiMatched->E());
248  // fPuppiCandidates->push_back(pCand);
249  } else {
250  pVec.SetPxPyPzE( 0, 0, 0, 0);
251  }
252  pCand.setP4(pVec);
253  puppiP4s.push_back( pVec );
254 
255  pCand.setSourceCandidatePtr( i0->sourceCandidatePtr(0) );
256  fPuppiCandidates->push_back(pCand);
257  }
258 
259  //Compute the modified p4s
260  edm::ValueMap<LorentzVector>::Filler p4PupFiller(*p4PupOut);
261  p4PupFiller.insert(hPFProduct,puppiP4s.begin(), puppiP4s.end() );
262  p4PupFiller.fill();
263 
264  iEvent.put(lPupOut);
265  iEvent.put(p4PupOut);
267  for(unsigned int ic=0, nc = oh->size(); ic < nc; ++ic) {
268  reco::CandidatePtr pkref( oh, ic );
269  values[ic] = pkref;
270 
271  }
272  std::auto_ptr<edm::ValueMap<reco::CandidatePtr> > pfMap_p(new edm::ValueMap<reco::CandidatePtr>());
274  filler.insert(hPFProduct, values.begin(), values.end());
275  filler.fill();
276  iEvent.put(pfMap_p);
277 
278 
281 
282  // all the different alphas per particle
283  // THE alpha per particle
284  std::auto_ptr<std::vector<double> > theAlphas(new std::vector<double>(fPuppiContainer->puppiAlphas()));
285  std::auto_ptr<std::vector<double> > theAlphasMed(new std::vector<double>(fPuppiContainer->puppiAlphasMed()));
286  std::auto_ptr<std::vector<double> > theAlphasRms(new std::vector<double>(fPuppiContainer->puppiAlphasRMS()));
287  std::auto_ptr<std::vector<double> > alphas(new std::vector<double>(fPuppiContainer->puppiRawAlphas()));
288  std::auto_ptr<double> nalgos(new double(fPuppiContainer->puppiNAlgos()));
289 
290  iEvent.put(alphas,"PuppiRawAlphas");
291  iEvent.put(nalgos,"PuppiNAlgos");
292  iEvent.put(theAlphas,"PuppiAlphas");
293  iEvent.put(theAlphasMed,"PuppiAlphasMed");
294  iEvent.put(theAlphasRms,"PuppiAlphasRms");
295  }
296 
297 }
298 
299 // ------------------------------------------------------------------------------------------
301 }
302 // ------------------------------------------------------------------------------------------
304 }
305 // ------------------------------------------------------------------------------------------
307 }
308 // ------------------------------------------------------------------------------------------
310 }
311 // ------------------------------------------------------------------------------------------
313 }
314 // ------------------------------------------------------------------------------------------
316 }
317 // ------------------------------------------------------------------------------------------
319  //The following says we do not know what parameters are allowed so do no validation
320  // Please change this to state exactly what you do use, even if it is no parameters
322  desc.setUnknown();
323  descriptions.addDefault(desc);
324 }
325 //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()
bool fUseWeightsNoLep
Definition: PuppiProducer.h:54
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
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:49
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:52
std::vector< RecoObj > fRecoObjCollection
Definition: PuppiProducer.h:58
virtual void setP4(const LorentzVector &p4)
set 4-momentum
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
std::vector< LorentzVector > LorentzVectorCollection
Definition: PuppiProducer.h:27
virtual void beginJob()
void setSourceCandidatePtr(const PFCandidatePtr &ptr)
Definition: PFCandidate.h:125
math::XYZTLorentzVector LorentzVector
Definition: PuppiProducer.h:26
const reco::VertexRef vertexRef() const
std::unique_ptr< PuppiContainer > fPuppiContainer
Definition: PuppiProducer.h:57
int charge
Definition: RecoObj.h:29
edm::EDGetTokenT< VertexCollection > tokenVertices_
Definition: PuppiProducer.h:45
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:115
const PVAssoc fromPV(size_t ipv=0) const
virtual double py() const
y coordinate of momentum vector
virtual double energy() const
energy
bool fPuppiForLeptons
Definition: PuppiProducer.h:50
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:53
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:85
reco::GsfTrackRef gsfTrackRef() const
Definition: PFCandidate.cc:471
edm::EDGetTokenT< CandidateView > tokenPFCandidates_
Definition: PuppiProducer.h:44
std::auto_ptr< PFOutputCollection > fPuppiCandidates
Definition: PuppiProducer.h:59
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
Definition: Run.h:41