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
21 //Main File
22 #include "fastjet/PseudoJet.hh"
24 
25 
26 
27 // ------------------------------------------------------------------------------------------
29  fUseDZ = iConfig.getParameter<bool>("UseDeltaZCut");
30  fDZCut = iConfig.getParameter<double>("DeltaZCut");
31  fPuppiContainer = std::unique_ptr<PuppiContainer> ( new PuppiContainer(iConfig) );
32 
34  = consumes<CandidateView>(iConfig.getParameter<edm::InputTag>("candName"));
36  = consumes<VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexName"));
37 
38 
39  produces<edm::ValueMap<float> > ("PuppiWeights");
40  produces<edm::ValueMap<LorentzVector> > ("PuppiP4s");
41  produces<PFOutputCollection>();
42 
43 
44 }
45 // ------------------------------------------------------------------------------------------
47 }
48 // ------------------------------------------------------------------------------------------
50 
51  // Get PFCandidate Collection
52  edm::Handle<CandidateView> hPFProduct;
53  iEvent.getByToken(tokenPFCandidates_,hPFProduct);
54  const CandidateView *pfCol = hPFProduct.product();
55 
56  // Get vertex collection w/PV as the first entry?
58  iEvent.getByToken(tokenVertices_,hVertexProduct);
59  const reco::VertexCollection *pvCol = hVertexProduct.product();
60 
61  //Fill the reco objects
62  fRecoObjCollection.clear();
63  for(CandidateView::const_iterator itPF = pfCol->begin(); itPF!=pfCol->end(); itPF++) {
64  RecoObj pReco;
65  pReco.pt = itPF->pt();
66  pReco.eta = itPF->eta();
67  pReco.phi = itPF->phi();
68  pReco.m = itPF->mass();
69  pReco.charge = itPF->charge();
70  const reco::Vertex *closestVtx = 0;
71  double pDZ = -9999;
72  double pD0 = -9999;
73  int pVtxId = -9999;
74  bool lFirst = true;
75  const pat::PackedCandidate *lPack = dynamic_cast<const pat::PackedCandidate*>(&(*itPF));
76  if(lPack == 0 ) {
77  const reco::PFCandidate *pPF = dynamic_cast<const reco::PFCandidate*>(&(*itPF));
78  for(reco::VertexCollection::const_iterator iV = pvCol->begin(); iV!=pvCol->end(); ++iV) {
79  if(lFirst) {
80  if ( pPF->trackRef().isNonnull() ) pDZ = pPF->trackRef() ->dz(iV->position());
81  else if ( pPF->gsfTrackRef().isNonnull() ) pDZ = pPF->gsfTrackRef()->dz(iV->position());
82  if ( pPF->trackRef().isNonnull() ) pD0 = pPF->trackRef() ->d0();
83  else if ( pPF->gsfTrackRef().isNonnull() ) pD0 = pPF->gsfTrackRef()->d0();
84  lFirst = false;
85  if(pDZ > -9999) pVtxId = 0;
86  }
87  if(iV->trackWeight(pPF->trackRef())>0) {
88  closestVtx = &(*iV);
89  break;
90  }
91  pVtxId++;
92  }
93  } else if(lPack->vertexRef().isNonnull() ) {
94  pDZ = lPack->dz();
95  pD0 = lPack->dxy();
96  closestVtx = &(*(lPack->vertexRef()));
97  pVtxId = (lPack->fromPV() != (pat::PackedCandidate::PVUsedInFit));
98  if( (lPack->fromPV() == pat::PackedCandidate::PVLoose) ||
99  (lPack->fromPV() == pat::PackedCandidate::PVTight) )
100  closestVtx = 0;
101  }
102  pReco.dZ = pDZ;
103  pReco.d0 = pD0;
104 
105  if(closestVtx == 0) pReco.vtxId = -1;
106  if(closestVtx != 0) pReco.vtxId = pVtxId;
107  //if(closestVtx != 0) pReco.vtxChi2 = closestVtx->trackWeight(itPF->trackRef());
108  //Set the id for Puppi Algo: 0 is neutral pfCandidate, id = 1 for particles coming from PV and id = 2 for charged particles from non-leading vertex
109  pReco.id = 0;
110 
111  if(closestVtx != 0 && pVtxId == 0 && fabs(pReco.charge) > 0) pReco.id = 1;
112  if(closestVtx != 0 && pVtxId > 0 && fabs(pReco.charge) > 0) pReco.id = 2;
113  //Add a dZ cut if wanted (this helps)
114  if(fUseDZ && pDZ > -9999 && closestVtx == 0 && (fabs(pDZ) < fDZCut) && fabs(pReco.charge) > 0) pReco.id = 1;
115  if(fUseDZ && pDZ > -9999 && closestVtx == 0 && (fabs(pDZ) > fDZCut) && fabs(pReco.charge) > 0) pReco.id = 2;
116 
117  //std::cout << "pVtxId = " << pVtxId << ", and charge = " << itPF->charge() << ", and closestVtx = " << closestVtx << ", and id = " << pReco.id << std::endl;
118 
119  fRecoObjCollection.push_back(pReco);
120  }
121  fPuppiContainer->initialize(fRecoObjCollection);
122 
123  //Compute the weights
124  const std::vector<double> lWeights = fPuppiContainer->puppiWeights();
125  //Fill it into the event
126  std::auto_ptr<edm::ValueMap<float> > lPupOut(new edm::ValueMap<float>());
127  edm::ValueMap<float>::Filler lPupFiller(*lPupOut);
128  lPupFiller.insert(hPFProduct,lWeights.begin(),lWeights.end());
129  lPupFiller.fill();
130 
131  // This is a dummy to access the "translate" method which is a
132  // non-static member function even though it doesn't need to be.
133  // Will fix in the future.
134  static const reco::PFCandidate dummySinceTranslateIsNotStatic;
135 
136  //Fill a new PF Candidate Collection and write out the ValueMap of the new p4s.
137  // Since the size of the ValueMap must be equal to the input collection, we need
138  // to search the "puppi" particles to find a match for each input. If none is found,
139  // the input is set to have a four-vector of 0,0,0,0
140  const std::vector<fastjet::PseudoJet> lCandidates = fPuppiContainer->puppiParticles();
142  std::auto_ptr<edm::ValueMap<LorentzVector> > p4PupOut(new edm::ValueMap<LorentzVector>());
143  LorentzVectorCollection puppiP4s;
144  for ( auto i0 = hPFProduct->begin(),
145  i0begin = hPFProduct->begin(),
146  i0end = hPFProduct->end(); i0 != i0end; ++i0 ) {
147  //for(unsigned int i0 = 0; i0 < lCandidates.size(); i0++) {
148  //reco::PFCandidate pCand;
149  auto id = dummySinceTranslateIsNotStatic.translatePdgIdToType(i0->pdgId());
150  reco::PFCandidate pCand( i0->charge(),
151  i0->p4(),
152  id );
153  LorentzVector pVec = i0->p4();
154  int val = i0 - i0begin;
155 
156  // Find the Puppi particle matched to the input collection using the "user_index" of the object.
157  auto puppiMatched = find_if( lCandidates.begin(), lCandidates.end(), [&val]( fastjet::PseudoJet const & i ){ return i.user_index() == val; } );
158  if ( puppiMatched != lCandidates.end() ) {
159  pVec.SetPxPyPzE(puppiMatched->px(),puppiMatched->py(),puppiMatched->pz(),puppiMatched->E());
160  } else {
161  pVec.SetPxPyPzE( 0, 0, 0, 0);
162  }
163  pCand.setP4(pVec);
164  puppiP4s.push_back( pVec );
165  fPuppiCandidates->push_back(pCand);
166  }
167 
168  //Compute the modified p4s
169  edm::ValueMap<LorentzVector>::Filler p4PupFiller(*p4PupOut);
170  p4PupFiller.insert(hPFProduct,puppiP4s.begin(), puppiP4s.end() );
171  p4PupFiller.fill();
172 
173  iEvent.put(lPupOut,"PuppiWeights");
174  iEvent.put(p4PupOut,"PuppiP4s");
175  iEvent.put(fPuppiCandidates);
176 }
177 
178 // ------------------------------------------------------------------------------------------
180 }
181 // ------------------------------------------------------------------------------------------
183 }
184 // ------------------------------------------------------------------------------------------
186 }
187 // ------------------------------------------------------------------------------------------
189 }
190 // ------------------------------------------------------------------------------------------
192 }
193 // ------------------------------------------------------------------------------------------
195 }
196 // ------------------------------------------------------------------------------------------
198  //The following says we do not know what parameters are allowed so do no validation
199  // Please change this to state exactly what you do use, even if it is no parameters
201  desc.setUnknown();
202  descriptions.addDefault(desc);
203 }
204 //define this as a plug-in
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:250
float d0
Definition: RecoObj.h:28
virtual void endJob()
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
std::vector< reco::PFCandidate > PFOutputCollection
Definition: PuppiProducer.h:31
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Definition: RecoObj.h:4
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:52
std::vector< RecoObj > fRecoObjCollection
Definition: PuppiProducer.h:52
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:51
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 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:113
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
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:44
std::auto_ptr< PFOutputCollection > fPuppiCandidates
Definition: PuppiProducer.h:53
const_iterator end() const
virtual void endRun(edm::Run &, edm::EventSetup const &)
virtual void beginLuminosityBlock(edm::LuminosityBlock &, edm::EventSetup const &)
virtual float dxy() const
dxy with respect to the PV ref
const PVAssoc fromPV() const
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
int vtxId
Definition: RecoObj.h:20
virtual float dz() const
dz with respect to the PV ref
Definition: Run.h:41