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