CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PATLostTracks.cc
Go to the documentation of this file.
1 #include <string>
2 
3 
26 
27 
28 namespace pat {
30  public:
31  explicit PATLostTracks(const edm::ParameterSet&);
33 
34  virtual void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
35 
36  private:
43  const double minPt_;
44  const double minHits_;
45  const double minPixelHits_;
46  };
47 }
48 
50  Cands_(consumes<reco::PFCandidateCollection>(iConfig.getParameter<edm::InputTag>("inputCandidates"))),
51  map_(consumes<edm::Association<pat::PackedCandidateCollection> >(iConfig.getParameter<edm::InputTag>("packedPFCandidates"))),
52  Tracks_(consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("inputTracks"))),
53  Vertices_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("secondaryVertices"))),
54  PV_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("primaryVertices"))),
55  PVOrigs_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("originalVertices"))),
56  minPt_(iConfig.getParameter<double>("minPt")),
57  minHits_(iConfig.getParameter<uint32_t>("minHits")),
58  minPixelHits_(iConfig.getParameter<uint32_t>("minPixelHits"))
59 {
60  produces< std::vector<reco::Track> > ();
61  produces< std::vector<pat::PackedCandidate> > ();
62  produces< edm::Association<pat::PackedCandidateCollection> > ();
63 }
64 
66 
68 
70  iEvent.getByToken( Cands_, cands );
71  std::vector<reco::Candidate>::const_iterator cand;
72 
74  iEvent.getByToken(map_,pf2pc);
75 
77  iEvent.getByToken( Tracks_, tracks );
78 
80  iEvent.getByToken( Vertices_, vertices );
81 
83  iEvent.getByToken( PV_, pvs );
84  reco::VertexRef PV(pvs.id());
85  if (!pvs->empty()) {
86  PV = reco::VertexRef(pvs, 0);
87  }
89  iEvent.getByToken( PVOrigs_, PVOrigs );
90  const reco::Vertex & PVOrig = (*PVOrigs)[0];
91 
92  std::auto_ptr< std::vector<reco::Track> > outPtrP( new std::vector<reco::Track> );
93  std::vector<int> used(tracks->size(),0);
94 
95  std::auto_ptr< std::vector<pat::PackedCandidate> > outPtrC( new std::vector<pat::PackedCandidate> );
96 
97  //Mark all tracks used in candidates
98  for(unsigned int ic=0, nc = cands->size(); ic < nc; ++ic) {
100  const reco::PFCandidate &cand=(*cands)[ic];
101  if (cand.charge()) {
102  if(cand.trackRef().isNonnull() && cand.trackRef().id() ==tracks.id() && (*pf2pc)[r]->numberOfHits() > 0) // also check if packed candidates are storing the tracks for this one
103  {
104  used[cand.trackRef().key()]=1;
105  }
106  }
107  }
108 
109  //Mark all tracks used in secondary vertices
110  for(unsigned int i=0; i < vertices->size(); i++){
111  const reco::Vertex & sv = (*vertices)[i];
112  for(reco::Vertex::trackRef_iterator it = sv.tracks_begin(),e=sv.tracks_end(); it!=e;it++){
113  if(used[it->key()]==0) used[it->key()]=2; // mark as white list
114  }
115  }
116 
117  std::vector<int> mapping(tracks->size(),-1);
118  int j=0;
119  for(unsigned int i=0; i < used.size(); i++)
120  {
121  const reco::Track & tr = (*tracks)[i];
122  if(used[i] == 2 ||
123  (used[i]==0 && tr.pt() > minPt_ && tr.numberOfValidHits() >= minHits_ && tr.hitPattern().numberOfValidPixelHits() >= minPixelHits_ )
124  )
125  {
126  outPtrP->push_back(tr);
127  reco::Candidate::PolarLorentzVector p4(tr.pt(),tr.eta(),tr.phi(),0.13957018);
128  outPtrC->push_back(pat::PackedCandidate(p4,tr.vertex(),tr.phi(),211*tr.charge(),PV));
129  outPtrC->back().setTrackProperties((*tracks)[i]);
130  if(PVOrig.trackWeight(edm::Ref<reco::TrackCollection>(tracks,i)) > 0.5) {
131  outPtrC->back().setAssociationQuality(pat::PackedCandidate::UsedInFitTight);
132  }
133 
134  mapping[i]=j;
135  j++;
136  }
137  }
138  iEvent.put(outPtrP);
140  std::auto_ptr<edm::Association<pat::PackedCandidateCollection> > tk2pc(new edm::Association<pat::PackedCandidateCollection>(oh ));
142  tk2pcFiller.insert(tracks, mapping.begin(), mapping.end());
143  tk2pcFiller.fill() ;
144  iEvent.put(tk2pc);
145 
146 }
147 
148 
149 using pat::PATLostTracks;
int i
Definition: DBlmapReader.cc:9
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
const double minPt_
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:44
ProductID id() const
Definition: HandleBase.cc:15
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< pat::PackedCandidate > PackedCandidateCollection
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:52
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:13
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:598
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
virtual void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
key_type key() const
Accessor for product key.
Definition: Ref.h:266
PATLostTracks(const edm::ParameterSet &)
const edm::EDGetTokenT< reco::VertexCollection > Vertices_
ProductID id() const
Accessor for product ID.
Definition: Ref.h:256
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:433
const double minHits_
const edm::EDGetTokenT< reco::VertexCollection > PVOrigs_
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead ...
Definition: TrackBase.h:640
int iEvent
Definition: GenABIO.cc:230
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:604
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
float trackWeight(const TrackBaseRef &r) const
returns the weight with which a Track has contributed to the vertex-fit.
double p4[4]
Definition: TauolaWrapper.h:92
double pt() const
track transverse momentum
Definition: TrackBase.h:574
virtual int charge() const
electric charge
Definition: LeafCandidate.h:91
int j
Definition: DBlmapReader.cc:9
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:773
edm::Ref< VertexCollection > VertexRef
persistent reference to a Vertex
Definition: VertexFwd.h:13
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
const edm::EDGetTokenT< edm::Association< pat::PackedCandidateCollection > > map_
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:411
tuple tracks
Definition: testEve_cfg.py:39
const edm::EDGetTokenT< reco::TrackCollection > Tracks_
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:37
int numberOfValidPixelHits() const
Definition: HitPattern.h:749
const double minPixelHits_
int charge() const
track electric charge
Definition: TrackBase.h:520
const edm::EDGetTokenT< reco::VertexCollection > PV_
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:39
const edm::EDGetTokenT< reco::PFCandidateCollection > Cands_
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:39