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 && std::abs(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 (std::abs(tmpdz) < curdz){
129  curdz = std::abs(tmpdz);
130  closestVtxForUnassociateds = pVtxId;
131  }
132  pVtxId++;
133 
134  }
135  int tmpFromPV = 0;
136  // mocking the miniAOD definitions
137  if (closestVtx != 0 && std::abs(pReco.charge) > 0 && pVtxId > 0) tmpFromPV = 0;
138  if (closestVtx != 0 && std::abs(pReco.charge) > 0 && pVtxId == 0) tmpFromPV = 3;
139  if (closestVtx == 0 && std::abs(pReco.charge) > 0 && closestVtxForUnassociateds == 0) tmpFromPV = 2;
140  if (closestVtx == 0 && std::abs(pReco.charge) > 0 && closestVtxForUnassociateds != 0) tmpFromPV = 1;
141  pReco.dZ = pDZ;
142  pReco.d0 = pD0;
143  pReco.id = 0;
144  if (std::abs(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 && (std::abs(pDZ) < fDZCut)) pReco.id = 1;
151  if (!fPuppiForLeptons && fUseDZ && (std::abs(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  pReco.dZ = pDZ;
161  pReco.d0 = pD0;
162 
163  pReco.id = 0;
164  if (std::abs(pReco.charge) == 0){ pReco.id = 0; }
165  if (std::abs(pReco.charge) > 0){
166  if (lPack->fromPV() == 0){ pReco.id = 2; } // 0 is associated to PU vertex
167  if (lPack->fromPV() == (pat::PackedCandidate::PVUsedInFit)){ pReco.id = 1; }
168  if (lPack->fromPV() == (pat::PackedCandidate::PVTight) || lPack->fromPV() == (pat::PackedCandidate::PVLoose)){
169  pReco.id = 0;
170  if (!fPuppiForLeptons && fUseDZ && (std::abs(pDZ) < fDZCut)) pReco.id = 1;
171  if (!fPuppiForLeptons && fUseDZ && (std::abs(pDZ) > fDZCut)) pReco.id = 2;
172  if (fPuppiForLeptons && lPack->fromPV() == (pat::PackedCandidate::PVLoose)) pReco.id = 2;
173  if (fPuppiForLeptons && lPack->fromPV() == (pat::PackedCandidate::PVTight)) pReco.id = 1;
174  }
175  }
176  }
177 
178  fRecoObjCollection.push_back(pReco);
179 
180  }
181 
182  fPuppiContainer->initialize(fRecoObjCollection);
183  fPuppiContainer->setNPV( npv );
184 
185  std::vector<double> lWeights;
186  std::vector<fastjet::PseudoJet> lCandidates;
187  if (!fUseExistingWeights){
188  //Compute the weights and get the particles
189  lWeights = fPuppiContainer->puppiWeights();
190  lCandidates = fPuppiContainer->puppiParticles();
191  }
192  else{
193  //Use the existing weights
194  int lPackCtr = 0;
195  for(CandidateView::const_iterator itPF = pfCol->begin(); itPF!=pfCol->end(); itPF++) {
196  const pat::PackedCandidate *lPack = dynamic_cast<const pat::PackedCandidate*>(&(*itPF));
197  float curpupweight = -1.;
198  if(lPack == 0 ) {
199  // throw error
200  throw edm::Exception(edm::errors::LogicError,"PuppiProducer: cannot get weights since inputs are not PackedCandidates");
201  }
202  else{
203  if (fUseWeightsNoLep){ curpupweight = lPack->puppiWeightNoLep(); }
204  else{ curpupweight = lPack->puppiWeight(); }
205  }
206  lWeights.push_back(curpupweight);
207  fastjet::PseudoJet curjet( curpupweight*lPack->px(), curpupweight*lPack->py(), curpupweight*lPack->pz(), curpupweight*lPack->energy());
208  curjet.set_user_index(lPackCtr);
209  lCandidates.push_back(curjet);
210  lPackCtr++;
211  }
212  }
213 
214  //Fill it into the event
215  std::unique_ptr<edm::ValueMap<float> > lPupOut(new edm::ValueMap<float>());
216  edm::ValueMap<float>::Filler lPupFiller(*lPupOut);
217  lPupFiller.insert(hPFProduct,lWeights.begin(),lWeights.end());
218  lPupFiller.fill();
219 
220  // This is a dummy to access the "translate" method which is a
221  // non-static member function even though it doesn't need to be.
222  // Will fix in the future.
223  static const reco::PFCandidate dummySinceTranslateIsNotStatic;
224 
225  // Fill a new PF/Packed Candidate Collection and write out the ValueMap of the new p4s.
226  // Since the size of the ValueMap must be equal to the input collection, we need
227  // to search the "puppi" particles to find a match for each input. If none is found,
228  // the input is set to have a four-vector of 0,0,0,0
231  std::unique_ptr<edm::ValueMap<LorentzVector> > p4PupOut(new edm::ValueMap<LorentzVector>());
232  LorentzVectorCollection puppiP4s;
233  std::vector<reco::CandidatePtr> values(hPFProduct->size());
234 
235  for ( auto i0 = hPFProduct->begin(),
236  i0begin = hPFProduct->begin(),
237  i0end = hPFProduct->end(); i0 != i0end; ++i0 ) {
238  std::unique_ptr<pat::PackedCandidate> pCand;
239  std::unique_ptr<reco::PFCandidate> pfCand;
241  const pat::PackedCandidate *cand = dynamic_cast<const pat::PackedCandidate*>(&(*i0));
242  if(!cand)
243  throw edm::Exception(edm::errors::LogicError,"PuppiProducer: inputs are not PackedCandidates");
244  pCand.reset( new pat::PackedCandidate(*cand) );
245  } else {
246  auto id = dummySinceTranslateIsNotStatic.translatePdgIdToType(i0->pdgId());
247  const reco::PFCandidate *cand = dynamic_cast<const reco::PFCandidate*>(&(*i0));
248  pfCand.reset( new reco::PFCandidate( cand ? *cand : reco::PFCandidate(i0->charge(), i0->p4(), id) ) );
249  }
250  LorentzVector pVec = i0->p4();
251  int val = i0 - i0begin;
252 
253  // Find the Puppi particle matched to the input collection using the "user_index" of the object.
254  auto puppiMatched = find_if( lCandidates.begin(), lCandidates.end(), [&val]( fastjet::PseudoJet const & i ){ return i.user_index() == val; } );
255  if ( puppiMatched != lCandidates.end() ) {
256  pVec.SetPxPyPzE(puppiMatched->px(),puppiMatched->py(),puppiMatched->pz(),puppiMatched->E());
257  } else {
258  pVec.SetPxPyPzE( 0, 0, 0, 0);
259  }
260  puppiP4s.push_back( pVec );
261 
263  pCand->setP4(pVec);
264  pCand->setSourceCandidatePtr( i0->sourceCandidatePtr(0) );
265  fPackedPuppiCandidates->push_back(*pCand);
266  } else {
267  pfCand->setP4(pVec);
268  pfCand->setSourceCandidatePtr( i0->sourceCandidatePtr(0) );
269  fPuppiCandidates->push_back(*pfCand);
270  }
271  }
272 
273  //Compute the modified p4s
274  edm::ValueMap<LorentzVector>::Filler p4PupFiller(*p4PupOut);
275  p4PupFiller.insert(hPFProduct,puppiP4s.begin(), puppiP4s.end() );
276  p4PupFiller.fill();
277 
278  iEvent.put(std::move(lPupOut));
279  iEvent.put(std::move(p4PupOut));
282  for(unsigned int ic=0, nc = oh->size(); ic < nc; ++ic) {
283  reco::CandidatePtr pkref( oh, ic );
284  values[ic] = pkref;
285  }
286  } else {
288  for(unsigned int ic=0, nc = oh->size(); ic < nc; ++ic) {
289  reco::CandidatePtr pkref( oh, ic );
290  values[ic] = pkref;
291  }
292  }
293  std::unique_ptr<edm::ValueMap<reco::CandidatePtr> > pfMap_p(new edm::ValueMap<reco::CandidatePtr>());
295  filler.insert(hPFProduct, values.begin(), values.end());
296  filler.fill();
297  iEvent.put(std::move(pfMap_p));
298 
299 
302 
303  // all the different alphas per particle
304  // THE alpha per particle
305  std::unique_ptr<std::vector<double> > theAlphas(new std::vector<double>(fPuppiContainer->puppiAlphas()));
306  std::unique_ptr<std::vector<double> > theAlphasMed(new std::vector<double>(fPuppiContainer->puppiAlphasMed()));
307  std::unique_ptr<std::vector<double> > theAlphasRms(new std::vector<double>(fPuppiContainer->puppiAlphasRMS()));
308  std::unique_ptr<std::vector<double> > alphas(new std::vector<double>(fPuppiContainer->puppiRawAlphas()));
309  std::unique_ptr<double> nalgos(new double(fPuppiContainer->puppiNAlgos()));
310 
311  iEvent.put(std::move(alphas),"PuppiRawAlphas");
312  iEvent.put(std::move(nalgos),"PuppiNAlgos");
313  iEvent.put(std::move(theAlphas),"PuppiAlphas");
314  iEvent.put(std::move(theAlphasMed),"PuppiAlphasMed");
315  iEvent.put(std::move(theAlphasRms),"PuppiAlphasRms");
316  }
317 
318 }
319 
320 // ------------------------------------------------------------------------------------------
322 }
323 // ------------------------------------------------------------------------------------------
325 }
326 // ------------------------------------------------------------------------------------------
328  //The following says we do not know what parameters are allowed so do no validation
329  // Please change this to state exactly what you do use, even if it is no parameters
331  desc.setUnknown();
332  descriptions.addDefault(desc);
333 }
334 //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 &)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
float puppiWeightNoLep() const
Weight from full PUPPI.
float d0
Definition: RecoObj.h:28
virtual void endJob()
bool fUseWeightsNoLep
Definition: PuppiProducer.h:50
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:45
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:55
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:54
int charge
Definition: RecoObj.h:29
edm::EDGetTokenT< VertexCollection > tokenVertices_
Definition: PuppiProducer.h:41
float dZ
Definition: RecoObj.h:27
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:438
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
const PVAssoc fromPV(size_t ipv=0) const
virtual double py() const
y coordinate of momentum vector
def move
Definition: eostools.py:510
std::unique_ptr< PackedOutputCollection > fPackedPuppiCandidates
Definition: PuppiProducer.h:57
virtual double energy() const
energy
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool fClonePackedCands
Definition: PuppiProducer.h:51
bool fPuppiForLeptons
Definition: PuppiProducer.h:46
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:229
float eta
Definition: RecoObj.h:17
bool fUseExistingWeights
Definition: PuppiProducer.h:49
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:476
edm::EDGetTokenT< CandidateView > tokenPFCandidates_
Definition: PuppiProducer.h:40
const_iterator end() const
virtual double px() const
x coordinate of momentum vector
virtual float dxy() const
dxy with respect to the PV ref
std::unique_ptr< PFOutputCollection > fPuppiCandidates
Definition: PuppiProducer.h:56
virtual double pz() const
z coordinate of momentum vector
virtual const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99