CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFCandidateMaker.cc
Go to the documentation of this file.
2 
4 
7 using namespace reco;
8 using namespace edm;
9 using namespace std;
10 
12 
13  PFCandidateCollection_ = iCollector.consumes<reco::PFCandidateCollection>(iConfig.getParameter<edm::InputTag>("pfCandsInputTag"));
14  PFElectrons_ = iCollector.consumes<edm::ValueMap<reco::PFCandidatePtr> >(iConfig.getParameter<edm::InputTag>("pfElectronsTag"));
15  TrackCollection_ = iCollector.consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("trackInputTag"));
16  thePVCollection_ = iCollector.consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("primaryVertexInputTag"));
17 
18 }
19 
20 void PFCandidateMaker::SetVars(HWW& hww, const edm::Event& iEvent, const edm::EventSetup& iSetup) {
21 
22  hww.Load_pfcands_p4();
23  hww.Load_pfcands_charge();
25  hww.Load_pfcands_vtxidx();
26  hww.Load_pfcands_trkidx();
28 
29  PFPileUpAlgo *pfPileUpAlgo_ = new PFPileUpAlgo();
30 
31  bool validToken;
32 
33  //get pfcandidates
35  Handle<PFCandidateCollection> pfCandidatesHandle;
36  validToken = iEvent.getByToken(PFCandidateCollection_, pfCandidatesHandle);
37  if(!validToken) return;
38  pfCandidates = pfCandidatesHandle.product();
39 
40  //get pfelectrons
42  Handle<PFCandMap> pfElectronsHandle;
43  validToken = iEvent.getByToken(PFElectrons_, pfElectronsHandle);
44  if(!validToken) return;
45  const PFCandMap *pfElectrons = pfElectronsHandle.product();
46 
47  // get tracks
49  validToken = iEvent.getByToken(TrackCollection_, track_h);
50  if(!validToken) return;
51 
52  // get vertices
54  validToken = iEvent.getByToken(thePVCollection_, vertex_h);
55  if(!validToken) return;
56  const reco::VertexCollection *vertices = vertex_h.product();
57 
58 
59  int iCand = 0;
60  for( PFCandidateCollection::const_iterator pf_it = pfCandidates->begin(); pf_it != pfCandidates->end(); pf_it++ ) {
61 
62  int pfflags = 0;
63  for( unsigned int i = 0; i < 17; i++ ) {
64  if(pf_it->flag((PFCandidate::Flags)i)) pfflags |= (1<<i);
65  }
66 
67  iCand++;
68 
69  hww.pfcands_p4() .push_back( LorentzVector(pf_it->px(), pf_it->py(), pf_it->pz(), pf_it->p()) );
70  hww.pfcands_charge() .push_back( pf_it->charge() );
71  hww.pfcands_particleId() .push_back( pf_it->translateTypeToPdgId(pf_it->particleId()) );
72 
73  //for charged pfcandidates, find corresponding track index
74  //here we take the track directly from PFCandidate::trackRef()
75  if( pf_it->charge() != 0 ){
76 
77  // match pf cand to a vertex
78  // using the no pileup algo
79  hww.pfcands_vtxidx().push_back(pfPileUpAlgo_->chargedHadronVertex(*vertices, *pf_it));
80 
81  reco::TrackRef pftrack = pf_it->trackRef();
82 
83  int trkidx = 0;
84  bool foundTrack = false;
85 
86  reco::TrackCollection::const_iterator tracks_end = track_h->end();
87 
88  for (reco::TrackCollection::const_iterator itrk = track_h->begin(); itrk != tracks_end; ++itrk) {
89 
90  reco::TrackRef trkref( track_h , itrk - track_h->begin() );
91 
92  if( pftrack.key() == trkref.key() ){
93 
94  //sanity check
95  float dpt = pftrack->pt() - trkref->pt();
96  if( fabs( dpt ) > 0.1 ){
97  edm::LogWarning("SanityCheck") << "Warning: pfcandidate track pt - matched track pt = " << dpt << ", possible mismatch";
98  }
99 
100  //found corresponding track
101  hww.pfcands_trkidx().push_back( trkidx );
102  foundTrack = true;
103  break;
104  }
105 
106  ++trkidx;
107  }
108 
109  if( !foundTrack ){
110  //no matched track found, set trkidx to -1
111  hww.pfcands_trkidx().push_back(-1);
112  }
113 
114  }else{
115  //neutral particle, set trkidx to -2
116  hww.pfcands_trkidx().push_back(-2);
117  // no vertex match
118  hww.pfcands_vtxidx().push_back(-2);
119  }
120 
121 
122  //find corresponding PFElectron index
123  if (pf_it->particleId() == PFCandidate::e){
124 
125  int index = -1;
126 
127  if (pf_it->gsfTrackRef().isNonnull()) {
128 
129  int pfGsfTkId = pf_it->gsfTrackRef().key();
130  unsigned int elsIndex = 0;
131  PFCandMap::const_iterator el_pit = pfElectrons->begin();
132  unsigned int nE = el_pit.size();
133  for(unsigned int iE = 0; iE < nE; ++iE){
134  const PFCandidatePtr& el_it = el_pit[iE];
135  if (el_it.isNull()) continue;
136  int elGsfTkId = -1;
137  if (el_it->gsfTrackRef().isNonnull()) elGsfTkId = el_it->gsfTrackRef().key();
138  if (elGsfTkId==pfGsfTkId) {
139  index = elsIndex;
140  break;
141 
142  }
143  elsIndex++;
144  }
145  }
146 
147  hww.pfcands_pfelsidx().push_back(index);
148 
149  } else {
150  hww.pfcands_pfelsidx().push_back(-2);
151  }
152 
153 
154  }//loop over candidate collection
155 
156 }
PFCandidateMaker(const edm::ParameterSet &, edm::ConsumesCollector)
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
void Load_pfcands_p4()
Definition: HWW.cc:1331
void Load_pfcands_pfelsidx()
Definition: HWW.cc:1340
key_type key() const
Definition: Ptr.h:169
void Load_pfcands_particleId()
Definition: HWW.cc:1337
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
std::vector< int > & pfcands_charge()
Definition: HWW.cc:745
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:13
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
std::vector< LorentzVector > & pfcands_p4()
Definition: HWW.cc:725
void Load_pfcands_trkidx()
Definition: HWW.cc:1334
key_type key() const
Accessor for product key.
Definition: Ref.h:266
std::pair< double, double > Point
Definition: CaloEllipse.h:18
std::vector< int > & pfcands_vtxidx()
Definition: HWW.cc:741
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
std::vector< int > & pfcands_trkidx()
Definition: HWW.cc:729
std::vector< int > & pfcands_pfelsidx()
Definition: HWW.cc:737
int iEvent
Definition: GenABIO.cc:230
std::vector< int > & pfcands_particleId()
Definition: HWW.cc:733
bool isNull() const
Checks for null.
Definition: Ptr.h:148
void Load_pfcands_charge()
Definition: HWW.cc:1346
Definition: HWW.h:12
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:152
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
T const * product() const
Definition: Handle.h:81
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
const_iterator begin() const
Definition: ValueMap.h:209
int chargedHadronVertex(const reco::VertexCollection &vertices, const reco::PFCandidate &pfcand) const
Definition: PFPileUpAlgo.cc:43
void Load_pfcands_vtxidx()
Definition: HWW.cc:1343
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< float > > XYZTLorentzVectorF
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:22
void SetVars(HWW &, const edm::Event &, const edm::EventSetup &)
edm::ValueMap< reco::PFCandidatePtr > PFCandMap
math::PtEtaPhiELorentzVectorF LorentzVector