CMS 3D CMS Logo

List of all members | Public Member Functions
ConversionHitChecker Class Reference

#include <ConversionHitChecker.h>

Public Member Functions

 ConversionHitChecker ()
 
std::pair< uint8_t, Measurement1DFloatnHitsBeforeVtx (const reco::TrackExtra &track, const reco::Vertex &vtx, float sigmaTolerance=3.0) const
 
uint8_t nSharedHits (const reco::Track &trk1, const reco::Track &trk2) const
 
 ~ConversionHitChecker ()
 

Detailed Description

Author
J.Bendavid
Version
Check hits along a Trajectory and count how many are before the vertex position. (taking into account the uncertainty in the vertex and hit positions. Also returns a the signed decay length and uncertainty from the closest hit on the track to the vertex position.

Definition at line 31 of file ConversionHitChecker.h.

Constructor & Destructor Documentation

ConversionHitChecker::ConversionHitChecker ( )
inline

Definition at line 35 of file ConversionHitChecker.h.

35 {}
ConversionHitChecker::~ConversionHitChecker ( )
inline

Member Function Documentation

std::pair< uint8_t, Measurement1DFloat > ConversionHitChecker::nHitsBeforeVtx ( const reco::TrackExtra track,
const reco::Vertex vtx,
float  sigmaTolerance = 3.0 
) const

Definition at line 9 of file ConversionHitChecker.cc.

References reco::Vertex::covariance(), listHistos::decayLength, dot(), h, hcalSimParameters_cfi::hb, mag2(), nhits, position, rpcPointValidation_cfi::recHit, reco::TrackExtraBase::recHitsBegin(), reco::TrackExtraBase::recHitsSize(), Measurement1DFloat::significance(), mathSSE::sqrt(), reco::TrackExtraBase::trajParams(), reco::Vertex::x(), reco::Vertex::y(), and reco::Vertex::z().

Referenced by ConversionProducer::buildCollection(), and ~ConversionHitChecker().

9  {
10 
11  // track hits are always inout
12 
13  GlobalPoint vtxPos(vtx.x(),vtx.y(),vtx.z());
14 
15  auto const & trajParams = track.trajParams();
16 
17  //iterate inside out, when distance to vertex starts increasing, we are at the closest hit
18  // the first (and last, btw) hit is always valid... (apparntly not..., conversion is different????)
19  auto hb = track.recHitsBegin();
20  unsigned int ih=0;
21  for(;ih<track.recHitsSize();++ih) {
22  auto recHit = *(hb+ih);
23  if (recHit->isValid()) break;
24  }
25  auto recHit = *(hb+ih);
26  unsigned int closest=ih;
27  auto globalPosition = recHit->surface()->toGlobal(trajParams[0].position());
28  auto distance2 = (vtxPos - globalPosition).mag2();
29  int nhits = 1;
30  for(unsigned int h=ih+1;h<track.recHitsSize();++h){
31 
32  //check if next valid hit is farther away from vertex than existing closest
33  auto nextHit = *(hb+h);
34  if (!nextHit->isValid() ) continue;
35  globalPosition = nextHit->surface()->toGlobal(trajParams[h].position());
36  auto nextDistance2 = (vtxPos - globalPosition).mag2();
37  if (nextDistance2 > distance2) break;
38 
39  distance2=nextDistance2;
40  ++nhits;
41  closest=h;
42  }
43 
44  //compute signed decaylength significance for closest hit and check if it is before the vertex
45  //if not then we need to subtract it from the count of hits before the vertex, since it has been implicitly included
46  recHit = *(hb+closest);
47  auto momDir = recHit->surface()->toGlobal(trajParams[closest].direction());
48  globalPosition = recHit->surface()->toGlobal(trajParams[closest].position());
49  float decayLengthHitToVtx = (vtxPos - globalPosition).dot(momDir);
50 
52  j[0] = momDir.x();
53  j[1] = momDir.y();
54  j[2] = momDir.z();
55  float vertexError2 = ROOT::Math::Similarity(j,vtx.covariance());
56  auto decayLenError = std::sqrt(vertexError2);
57 
58  Measurement1DFloat decayLength(decayLengthHitToVtx,decayLenError);
59 
60  if (decayLength.significance() < sigmaTolerance) { //decay length is not (significantly) positive, so hit is consistent with the vertex position or late
61  //subtract it from wrong hits count
62  --nhits;
63  }
64 
65  return std::pair<unsigned int,Measurement1DFloat>(nhits,decayLength);
66 
67 }
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
double y() const
y coordinate
Definition: Vertex.h:113
double covariance(int i, int j) const
(i, j)-th element of error matrix, i, j = 0, ... 2
Definition: Vertex.h:130
unsigned int recHitsSize() const
number of RecHits
T sqrt(T t)
Definition: SSEVec.h:18
TrajParams const & trajParams() const
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
ROOT::Math::SVector< double, 3 > AlgebraicVector3
double z() const
z coordinate
Definition: Vertex.h:115
trackingRecHit_iterator recHitsBegin() const
first iterator over RecHits
double x() const
x coordinate
Definition: Vertex.h:111
T dot(const Basic3DVector &v) const
Scalar product, or "dot" product, with a vector of same type.
static int position[264][3]
Definition: ReadPGInfo.cc:509
uint8_t ConversionHitChecker::nSharedHits ( const reco::Track trk1,
const reco::Track trk2 
) const

Definition at line 69 of file ConversionHitChecker.cc.

References TrackingRecHit::isValid(), reco::Track::recHitsBegin(), reco::Track::recHitsEnd(), TrackingRecHit::sharesInput(), and TrackingRecHit::some.

Referenced by ConversionProducer::buildCollection(), and ~ConversionHitChecker().

69  {
70 
71  uint8_t nShared = 0;
72 
73  for (trackingRecHit_iterator iHit1 = trk1.recHitsBegin(); iHit1 != trk1.recHitsEnd(); ++iHit1) {
74  const TrackingRecHit *hit1 = (*iHit1);
75  if (hit1->isValid()) {
76  for (trackingRecHit_iterator iHit2 = trk2.recHitsBegin(); iHit2 != trk2.recHitsEnd(); ++iHit2) {
77  const TrackingRecHit *hit2 = (*iHit2);
78  if (hit2->isValid() && hit1->sharesInput(hit2,TrackingRecHit::some)) {
79  ++nShared;
80  }
81  }
82  }
83  }
84 
85  return nShared;
86 
87 }
virtual bool sharesInput(const TrackingRecHit *other, SharedInputType what) const
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:106
bool isValid() const
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:111