CMS 3D CMS Logo

PATLostTracks.cc
Go to the documentation of this file.
14 
15 namespace {
16  bool passesQuality(const reco::Track& trk, const std::vector<reco::TrackBase::TrackQuality>& allowedQuals) {
17  for (const auto& qual : allowedQuals) {
18  if (trk.quality(qual))
19  return true;
20  }
21  return false;
22  }
23 } // namespace
24 
25 namespace pat {
27  public:
28  explicit PATLostTracks(const edm::ParameterSet&);
29  ~PATLostTracks() override;
30 
31  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
32 
33  private:
35  bool passTrkCuts(const reco::Track& tr) const;
36  void addPackedCandidate(std::vector<pat::PackedCandidate>& cands,
37  const reco::TrackRef& trk,
38  const reco::VertexRef& pvSlimmed,
39  const reco::VertexRefProd& pvSlimmedColl,
40  const TrkStatus& trkStatus,
41  const pat::PackedCandidate::PVAssociationQuality& pvAssocQuality,
43  std::pair<int, pat::PackedCandidate::PVAssociationQuality> associateTrkToVtx(const reco::VertexCollection& vertices,
44  const reco::TrackRef& trk) const;
45 
46  private:
55  const double minPt_;
56  const double minHits_;
57  const double minPixelHits_;
58  const double minPtToStoreProps_;
60  const int covarianceVersion_;
61  const std::vector<int> covariancePackingSchemas_;
62  std::vector<reco::TrackBase::TrackQuality> qualsToAutoAccept_;
70  const bool useLegacySetup_;
71  const bool xiSelection_;
72  const double xiMassCut_;
73  };
74 } // namespace pat
75 
77  : cands_(consumes<reco::PFCandidateCollection>(iConfig.getParameter<edm::InputTag>("inputCandidates"))),
78  map_(consumes<edm::Association<pat::PackedCandidateCollection>>(
79  iConfig.getParameter<edm::InputTag>("packedPFCandidates"))),
80  tracks_(consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("inputTracks"))),
81  vertices_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("secondaryVertices"))),
82  kshorts_(consumes<reco::VertexCompositeCandidateCollection>(iConfig.getParameter<edm::InputTag>("kshorts"))),
83  lambdas_(consumes<reco::VertexCompositeCandidateCollection>(iConfig.getParameter<edm::InputTag>("lambdas"))),
84  pv_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("primaryVertices"))),
85  pvOrigs_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("originalVertices"))),
86  minPt_(iConfig.getParameter<double>("minPt")),
87  minHits_(iConfig.getParameter<uint32_t>("minHits")),
88  minPixelHits_(iConfig.getParameter<uint32_t>("minPixelHits")),
89  minPtToStoreProps_(iConfig.getParameter<double>("minPtToStoreProps")),
90  minPtToStoreLowQualityProps_(iConfig.getParameter<double>("minPtToStoreLowQualityProps")),
91  covarianceVersion_(iConfig.getParameter<int>("covarianceVersion")),
92  covariancePackingSchemas_(iConfig.getParameter<std::vector<int>>("covariancePackingSchemas")),
93  muons_(consumes<reco::MuonCollection>(iConfig.getParameter<edm::InputTag>("muons"))),
94  passThroughCut_(iConfig.getParameter<std::string>("passThroughCut")),
95  maxDzForPrimaryAssignment_(
96  iConfig.getParameter<edm::ParameterSet>("pvAssignment").getParameter<double>("maxDzForPrimaryAssignment")),
97  maxDzSigForPrimaryAssignment_(
98  iConfig.getParameter<edm::ParameterSet>("pvAssignment").getParameter<double>("maxDzSigForPrimaryAssignment")),
99  maxDzErrorForPrimaryAssignment_(iConfig.getParameter<edm::ParameterSet>("pvAssignment")
100  .getParameter<double>("maxDzErrorForPrimaryAssignment")),
101  maxDxyForNotReconstructedPrimary_(iConfig.getParameter<edm::ParameterSet>("pvAssignment")
102  .getParameter<double>("maxDxyForNotReconstructedPrimary")),
103  maxDxySigForNotReconstructedPrimary_(iConfig.getParameter<edm::ParameterSet>("pvAssignment")
104  .getParameter<double>("maxDxySigForNotReconstructedPrimary")),
105  useLegacySetup_(iConfig.getParameter<bool>("useLegacySetup")),
106  xiSelection_(iConfig.getParameter<bool>("xiSelection")),
107  xiMassCut_(iConfig.getParameter<double>("xiMassCut")) {
108  std::vector<std::string> trkQuals(iConfig.getParameter<std::vector<std::string>>("qualsToAutoAccept"));
110  trkQuals.begin(), trkQuals.end(), std::back_inserter(qualsToAutoAccept_), reco::TrackBase::qualityByName);
111 
113  qualsToAutoAccept_.end()) {
114  std::ostringstream msg;
115  msg << " PATLostTracks has a quality requirement which resolves to undefQuality. This usually means a typo and is "
116  "therefore treated a config error\nquality requirements:\n ";
117  for (const auto& trkQual : trkQuals)
118  msg << trkQual << " ";
119  throw cms::Exception("Configuration") << msg.str();
120  }
121 
122  produces<std::vector<reco::Track>>();
123  produces<std::vector<pat::PackedCandidate>>();
124  produces<std::vector<pat::PackedCandidate>>("eleTracks");
125  produces<edm::Association<pat::PackedCandidateCollection>>();
126 }
127 
129 
132  iEvent.getByToken(cands_, cands);
133 
135  iEvent.getByToken(map_, pf2pc);
136 
138  iEvent.getByToken(tracks_, tracks);
139 
141  iEvent.getByToken(vertices_, vertices);
142 
144  iEvent.getByToken(muons_, muons);
145 
147  iEvent.getByToken(kshorts_, kshorts);
149  iEvent.getByToken(lambdas_, lambdas);
150 
152  iEvent.getByToken(pv_, pvs);
153  reco::VertexRef pv(pvs.id());
154  reco::VertexRefProd pvRefProd(pvs);
156  iEvent.getByToken(pvOrigs_, pvOrigs);
157 
158  auto outPtrTrks = std::make_unique<std::vector<reco::Track>>();
159  auto outPtrTrksAsCands = std::make_unique<std::vector<pat::PackedCandidate>>();
160  auto outPtrEleTrksAsCands = std::make_unique<std::vector<pat::PackedCandidate>>();
161 
162  std::vector<TrkStatus> trkStatus(tracks->size(), TrkStatus::NOTUSED);
163  //Mark all tracks used in candidates
164  //check if packed candidates are storing the tracks by seeing if number of hits >0
165  //currently we dont use that information though
166  //electrons will never store their track (they store the GSF track)
167  for (unsigned int ic = 0, nc = cands->size(); ic < nc; ++ic) {
169  const reco::PFCandidate& cand = (*cands)[ic];
170  if (cand.charge() && cand.trackRef().isNonnull() && cand.trackRef().id() == tracks.id()) {
171  if (cand.pdgId() == 11)
172  trkStatus[cand.trackRef().key()] = TrkStatus::PFELECTRON;
173  else if (cand.pdgId() == -11)
174  trkStatus[cand.trackRef().key()] = TrkStatus::PFPOSITRON;
175  else if ((*pf2pc)[r]->numberOfHits() > 0)
176  trkStatus[cand.trackRef().key()] = TrkStatus::PFCAND;
177  else
178  trkStatus[cand.trackRef().key()] = TrkStatus::PFCANDNOTRKPROPS;
179  }
180  }
181 
182  //Mark all tracks used in secondary vertices
183  for (const auto& secVert : *vertices) {
184  for (auto trkIt = secVert.tracks_begin(); trkIt != secVert.tracks_end(); trkIt++) {
185  if (trkStatus[trkIt->key()] == TrkStatus::NOTUSED)
186  trkStatus[trkIt->key()] = TrkStatus::VTX;
187  }
188  }
189  for (const auto& v0 : *kshorts) {
190  for (size_t dIdx = 0; dIdx < v0.numberOfDaughters(); dIdx++) {
191  size_t key = (dynamic_cast<const reco::RecoChargedCandidate*>(v0.daughter(dIdx)))->track().key();
192  if (trkStatus[key] == TrkStatus::NOTUSED)
193  trkStatus[key] = TrkStatus::VTX;
194  }
195  }
196  for (const auto& v0 : *lambdas) {
197  double protonCharge = 0;
198  for (size_t dIdx = 0; dIdx < v0.numberOfDaughters(); dIdx++) {
199  size_t key = (dynamic_cast<const reco::RecoChargedCandidate*>(v0.daughter(dIdx)))->track().key();
200  if (trkStatus[key] == TrkStatus::NOTUSED)
201  trkStatus[key] = TrkStatus::VTX;
202  protonCharge += v0.daughter(dIdx)->charge() * v0.daughter(dIdx)->momentum().mag2();
203  }
204  if (xiSelection_) {
205  // selecting potential Xi- -> Lambda pi candidates
206  math::XYZTLorentzVector p4Lambda(v0.px(), v0.py(), v0.pz(), sqrt(v0.momentum().mag2() + v0.mass() * v0.mass()));
207 
208  for (unsigned int trkIndx = 0; trkIndx < tracks->size(); trkIndx++) {
209  reco::TrackRef trk(tracks, trkIndx);
210  if ((*trk).charge() * protonCharge < 0 || trkStatus[trkIndx] != TrkStatus::NOTUSED)
211  continue;
212 
214  trk->px(), trk->py(), trk->pz(), sqrt(trk->momentum().mag2() + 0.01947995518)); // pion mass ^2
215  if ((p4Lambda + p4pi).M() < xiMassCut_)
216  trkStatus[trkIndx] = TrkStatus::VTX; // selecting potential Xi- candidates
217  }
218  }
219  }
220  std::vector<int> mapping(tracks->size(), -1);
221  int lostTrkIndx = 0;
222  for (unsigned int trkIndx = 0; trkIndx < tracks->size(); trkIndx++) {
223  reco::TrackRef trk(tracks, trkIndx);
224  if (trkStatus[trkIndx] == TrkStatus::VTX || (trkStatus[trkIndx] == TrkStatus::NOTUSED && passTrkCuts(*trk))) {
225  outPtrTrks->emplace_back(*trk);
226  //association to PV
227  std::pair<int, pat::PackedCandidate::PVAssociationQuality> pvAsso = associateTrkToVtx(*pvOrigs, trk);
228  const reco::VertexRef& pvOrigRef = reco::VertexRef(pvOrigs, pvAsso.first);
229  if (pvOrigRef.isNonnull()) {
230  pv = reco::VertexRef(pvs, pvOrigRef.key()); // WARNING: assume the PV slimmer is keeping same order
231  } else if (!pvs->empty()) {
232  pv = reco::VertexRef(pvs, 0);
233  }
234  addPackedCandidate(*outPtrTrksAsCands, trk, pv, pvRefProd, trkStatus[trkIndx], pvAsso.second, muons);
235 
236  //for creating the reco::Track -> pat::PackedCandidate map
237  //not done for the lostTrack:eleTracks collection
238  mapping[trkIndx] = lostTrkIndx;
239  lostTrkIndx++;
240  } else if ((trkStatus[trkIndx] == TrkStatus::PFELECTRON || trkStatus[trkIndx] == TrkStatus::PFPOSITRON) &&
241  passTrkCuts(*trk)) {
242  //association to PV
243  std::pair<int, pat::PackedCandidate::PVAssociationQuality> pvAsso = associateTrkToVtx(*pvOrigs, trk);
244  const reco::VertexRef& pvOrigRef = reco::VertexRef(pvOrigs, pvAsso.first);
245  if (pvOrigRef.isNonnull()) {
246  pv = reco::VertexRef(pvs, pvOrigRef.key()); // WARNING: assume the PV slimmer is keeping same order
247  } else if (!pvs->empty()) {
248  pv = reco::VertexRef(pvs, 0);
249  }
250  addPackedCandidate(*outPtrEleTrksAsCands, trk, pv, pvRefProd, trkStatus[trkIndx], pvAsso.second, muons);
251  }
252  }
253 
254  iEvent.put(std::move(outPtrTrks));
255  iEvent.put(std::move(outPtrEleTrksAsCands), "eleTracks");
257  auto tk2pc = std::make_unique<edm::Association<pat::PackedCandidateCollection>>(oh);
259  tk2pcFiller.insert(tracks, mapping.begin(), mapping.end());
260  tk2pcFiller.fill();
261  iEvent.put(std::move(tk2pc));
262 }
263 
265  const bool passTrkHits = tr.pt() > minPt_ && tr.numberOfValidHits() >= minHits_ &&
266  tr.hitPattern().numberOfValidPixelHits() >= minPixelHits_;
267  const bool passTrkQual = passesQuality(tr, qualsToAutoAccept_);
268 
269  return passTrkHits || passTrkQual || passThroughCut_(tr);
270 }
271 
272 void pat::PATLostTracks::addPackedCandidate(std::vector<pat::PackedCandidate>& cands,
273  const reco::TrackRef& trk,
274  const reco::VertexRef& pvSlimmed,
275  const reco::VertexRefProd& pvSlimmedColl,
276  const pat::PATLostTracks::TrkStatus& trkStatus,
277  const pat::PackedCandidate::PVAssociationQuality& pvAssocQuality,
279  const float mass = 0.13957018;
280 
281  int id = 211 * trk->charge();
282  if (trkStatus == TrkStatus::PFELECTRON)
283  id = 11;
284  else if (trkStatus == TrkStatus::PFPOSITRON)
285  id = -11;
286 
287  // assign the proper pdgId for tracks that are reconstructed as a muon
288  const reco::Muon* muon(nullptr);
289  for (auto& mu : *muons) {
290  if (reco::TrackRef(mu.innerTrack()) == trk) {
291  id = -13 * trk->charge();
292  muon = &mu;
293  break;
294  }
295  }
296 
298  int nlost = trk->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS);
299  if (nlost == 0) {
300  if (trk->hitPattern().hasValidHitInPixelLayer(PixelSubdetector::SubDetector::PixelBarrel, 1)) {
302  }
303  } else {
305  }
306 
307  reco::Candidate::PolarLorentzVector p4(trk->pt(), trk->eta(), trk->phi(), mass);
308  cands.emplace_back(
309  pat::PackedCandidate(p4, trk->vertex(), trk->pt(), trk->eta(), trk->phi(), id, pvSlimmedColl, pvSlimmed.key()));
310 
311  cands.back().setTrackHighPurity(trk->quality(reco::TrackBase::highPurity));
312  cands.back().setCovarianceVersion(covarianceVersion_);
313  cands.back().setLostInnerHits(lostHits);
314  if (trk->pt() > minPtToStoreProps_ || trkStatus == TrkStatus::VTX) {
315  if (useLegacySetup_ || std::abs(id) == 11 || trkStatus == TrkStatus::VTX) {
316  cands.back().setTrackProperties(*trk, covariancePackingSchemas_[4], covarianceVersion_);
317  } else {
318  if (trk->hitPattern().numberOfValidPixelHits() > 0) {
319  cands.back().setTrackProperties(
320  *trk, covariancePackingSchemas_[0], covarianceVersion_); // high quality with pixels
321  } else {
322  cands.back().setTrackProperties(
323  *trk, covariancePackingSchemas_[1], covarianceVersion_); // high quality without pixels
324  }
325  }
326  } else if (!useLegacySetup_ && trk->pt() > minPtToStoreLowQualityProps_) {
327  if (trk->hitPattern().numberOfValidPixelHits() > 0) {
328  cands.back().setTrackProperties(
329  *trk, covariancePackingSchemas_[2], covarianceVersion_); // low quality with pixels
330  } else {
331  cands.back().setTrackProperties(
332  *trk, covariancePackingSchemas_[3], covarianceVersion_); // low quality without pixels
333  }
334  }
335  cands.back().setAssociationQuality(pvAssocQuality);
336 
337  if (muon)
338  cands.back().setMuonID(muon->isStandAloneMuon(), muon->isGlobalMuon());
339 }
340 
341 std::pair<int, pat::PackedCandidate::PVAssociationQuality> pat::PATLostTracks::associateTrkToVtx(
342  const reco::VertexCollection& vertices, const reco::TrackRef& trk) const {
343  //For legacy setup check only if the track is used in fit of the PV, i.e. vertices[0],
344  //and associate quality if weight > 0.5. Otherwise return invalid vertex index (-1)
345  //and default quality flag (NotReconstructedPrimary = 0)
346  if (useLegacySetup_) {
347  float w = vertices[0].trackWeight(trk);
348  if (w > 0.5) {
349  return std::pair<int, pat::PackedCandidate::PVAssociationQuality>(0, pat::PackedCandidate::UsedInFitTight);
350  } else if (w > 0.) {
351  return std::pair<int, pat::PackedCandidate::PVAssociationQuality>(0,
353  } else {
354  return std::pair<int, pat::PackedCandidate::PVAssociationQuality>(-1,
356  }
357  }
358 
359  //Inspired by CommonTools/RecoAlgos/interface/PrimaryVertexAssignment.h
360  //but without specific association for secondaries in jets and option to use timing
361 
362  int iVtxMaxWeight = -1;
363  int iVtxMinDzDist = -1;
364  size_t idx = 0;
365  float maxWeight = 0;
367  double minDzSig = std::numeric_limits<double>::max();
368  for (auto const& vtx : vertices) {
369  float w = vtx.trackWeight(trk);
370  double dz = std::abs(trk->dz(vtx.position()));
371  double dzSig = dz / trk->dzError();
372  if (w > maxWeight) {
373  maxWeight = w;
374  iVtxMaxWeight = idx;
375  }
376  if (dzSig < minDzSig) {
377  minDzSig = dzSig;
378  minDz = dz;
379  iVtxMinDzDist = idx;
380  }
381  idx++;
382  }
383  // vertex in which fit the track was used
384  if (iVtxMaxWeight >= 0) {
385  if (maxWeight > 0.5) {
386  return std::pair<int, pat::PackedCandidate::PVAssociationQuality>(iVtxMaxWeight,
388  } else {
389  return std::pair<int, pat::PackedCandidate::PVAssociationQuality>(iVtxMaxWeight,
391  }
392  }
393  // vertex "closest in Z" with tight cuts (targetting primary particles)
394  if (minDz < maxDzForPrimaryAssignment_) {
395  const double add_cov = vertices[iVtxMinDzDist].covariance(2, 2);
396  const double dzErr = sqrt(trk->dzError() * trk->dzError() + add_cov);
397  if (minDz / dzErr < maxDzSigForPrimaryAssignment_ && trk->dzError() < maxDzErrorForPrimaryAssignment_) {
398  return std::pair<int, pat::PackedCandidate::PVAssociationQuality>(iVtxMinDzDist,
400  }
401  }
402  // if the track is not compatible with other PVs but is compatible with the BeamSpot, we may simply have not reco'ed the PV!
403  // we still point it to the closest in Z, but flag it as possible orphan-primary
404  if (!vertices.empty() && std::abs(trk->dxy(vertices[0].position())) < maxDxyForNotReconstructedPrimary_ &&
405  std::abs(trk->dxy(vertices[0].position()) / trk->dxyError()) < maxDxySigForNotReconstructedPrimary_)
406  return std::pair<int, pat::PackedCandidate::PVAssociationQuality>(iVtxMinDzDist,
408  // for tracks not associated to any PV return the closest in dz
409  return std::pair<int, pat::PackedCandidate::PVAssociationQuality>(iVtxMinDzDist, pat::PackedCandidate::OtherDeltaZ);
410 }
411 
412 using pat::PATLostTracks;
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:552
int numberOfValidPixelHits() const
Definition: HitPattern.h:831
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::pair< int, pat::PackedCandidate::PVAssociationQuality > associateTrkToVtx(const reco::VertexCollection &vertices, const reco::TrackRef &trk) const
std::vector< VertexCompositeCandidate > VertexCompositeCandidateCollection
collection of Candidate objects
const double minPt_
const double minPtToStoreProps_
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:798
T w() const
~PATLostTracks() override
std::vector< l1t::PFCandidate > PFCandidateCollection
Definition: PFCandidate.h:82
std::vector< pat::PackedCandidate > PackedCandidateCollection
const double maxDzSigForPrimaryAssignment_
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
const double minPtToStoreLowQualityProps_
PATLostTracks(const edm::ParameterSet &)
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
const double maxDxyForNotReconstructedPrimary_
std::vector< Vertex > VertexCollection
Definition: Vertex.h:31
const int covarianceVersion_
key_type key() const
Accessor for product key.
Definition: Ref.h:250
Definition: HeavyIon.h:7
const double minHits_
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
double pt() const
track transverse momentum
Definition: TrackBase.h:637
const edm::EDGetTokenT< reco::PFCandidateCollection > cands_
int iEvent
Definition: GenABIO.cc:224
const double maxDzErrorForPrimaryAssignment_
const edm::EDGetTokenT< reco::MuonCollection > muons_
const bool useLegacySetup_
T sqrt(T t)
Definition: SSEVec.h:19
const edm::EDGetTokenT< reco::VertexCompositeCandidateCollection > lambdas_
def pv(vc)
Definition: MetAnalyzer.py:7
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const double xiMassCut_
edm::Ref< VertexCollection > VertexRef
persistent reference to a Vertex
Definition: VertexFwd.h:13
const double maxDzForPrimaryAssignment_
const double maxDxySigForNotReconstructedPrimary_
LostInnerHits
Enumerator specifying the.
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:126
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:504
const edm::EDGetTokenT< reco::TrackCollection > tracks_
const edm::EDGetTokenT< reco::VertexCollection > vertices_
const std::vector< int > covariancePackingSchemas_
const bool xiSelection_
tuple msg
Definition: mps_check.py:286
void addPackedCandidate(std::vector< pat::PackedCandidate > &cands, const reco::TrackRef &trk, const reco::VertexRef &pvSlimmed, const reco::VertexRefProd &pvSlimmedColl, const TrkStatus &trkStatus, const pat::PackedCandidate::PVAssociationQuality &pvAssocQuality, edm::Handle< reco::MuonCollection > muons) const
const edm::EDGetTokenT< reco::VertexCompositeCandidateCollection > kshorts_
std::vector< Muon > MuonCollection
Definition: Muon.h:35
std::vector< reco::TrackBase::TrackQuality > qualsToAutoAccept_
StringCutObjectSelector< reco::Track, false > passThroughCut_
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
fixed size matrix
HLT enums.
bool passTrkCuts(const reco::Track &tr) const
const double minPixelHits_
const edm::EDGetTokenT< reco::VertexCollection > pv_
def move(src, dest)
Definition: eostools.py:511
const edm::EDGetTokenT< reco::VertexCollection > pvOrigs_
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:38
unsigned transform(const HcalDetId &id, unsigned transformCode)