CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GetTrackTrajInfo.cc
Go to the documentation of this file.
2 
6 
7 // To convert detId to subdet/layer number.
10 //#include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
13 
16 
20 
23 
25 
27 
28 std::vector< GetTrackTrajInfo::Result > GetTrackTrajInfo::analyze(const edm::EventSetup& iSetup, const reco::Track& track)
29 {
30  // Determine the track trajectory and detLayer at each layer that the track produces a hit in.
31 
32  std::vector< GetTrackTrajInfo::Result > results;
33 
34  // Initialise Tracker geometry info (not sufficient to do this only on first call).
36  iSetup.get<TrackerRecoGeometryRecord>().get( tracker );
37 
38  // This is also needed to extrapolate amongst the tracker layers.
40  iSetup.get<NavigationSchoolRecord>().get("SimpleNavigationSchool",theSchool);
41  NavigationSetter junk(*theSchool);
42 
43  // Get the magnetic field and use it to define a propagator for extrapolating the track trajectory.
45  iSetup.get<IdealMagneticFieldRecord>().get(magField);
47 
48  // This is used to check if a track is compatible with crossing a sensor.
49  // Use +3.0 rather than default -3.0 here, so hit defined as inside acceptance if
50  // no more than 3*sigma outside detector edge, as opposed to more than 3*sigma inside detector edge.
51  Chi2MeasurementEstimator estimator(30.,3.0);
52 
53  // Convert track to transientTrack, and hence to TSOS at its point of closest approach to beam.
54  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",trkTool_); // Needed for vertex fits
55  reco::TransientTrack t_trk = trkTool_->build(track);
56  TrajectoryStateOnSurface initTSOS = t_trk.impactPointState();
57  LogDebug("GTTI")<<"TRACK TSOS POS: x="<<initTSOS.globalPosition().x()<<" y="<<initTSOS.globalPosition().y()<<" z="<<initTSOS.globalPosition().z();
58 
59  // Note if the track is going into +ve or -ve z.
60  // This is only used to guess if the track is more likely to have hit a +ve rather than a -ve endcap
61  // disk. Since the +ve and -ve disks are a long way apart, this approximate method is good enough.
62  // More precise would be to check both possiblities and see which one (if either) the track crosses
63  // using detLayer::compatible().
64  bool posSide = track.eta() > 0;
65 
66  // Get hit patterns of this track
67  const reco::HitPattern& hp = track.hitPattern();
68 
69  // Loop over info for each hit
70  // N.B. Hits are sorted according to increasing distance from origin by
71  // RecoTracker/TrackProducer/src/TrackProducerBase.cc
72  for (int i = 0; i < hp.numberOfHits(); i++) {
73  uint32_t hit = hp.getHitPattern(i);
74  if (hp.trackerHitFilter(hit) && hp.validHitFilter(hit)) {
75  uint32_t subDet = hp.getSubStructure(hit);
76  uint32_t layer = hp.getLayer(hit);
77  // subdet: PixelBarrel=1, PixelEndcap=2, TIB=3, TID=4, TOB=5, TEC=6
78  LogDebug("GTTI")<<" hit in subdet="<<subDet<<" layer="<<layer;
79 
80  // Get corresponding DetLayer object (based on code in GeometricSearchTracker::idToLayer(...)
81  const DetLayer* detLayer = 0;
82  if (subDet == StripSubdetector::TIB) {
83  detLayer = tracker->tibLayers()[layer - 1];
84  } else if (subDet == StripSubdetector::TOB) {
85  detLayer = tracker->tobLayers()[layer - 1];
86  } else if (subDet == StripSubdetector::TID) {
87  detLayer = posSide ? tracker->posTidLayers()[layer - 1] : tracker->negTidLayers()[layer - 1];
88  } else if (subDet == StripSubdetector::TEC) {
89  detLayer = posSide ? tracker->posTecLayers()[layer - 1] : tracker->negTecLayers()[layer - 1];
90  } else if (subDet == PixelSubdetector::PixelBarrel) {
91  detLayer = tracker->pixelBarrelLayers()[layer - 1];
92  } else if (subDet == PixelSubdetector::PixelEndcap) {
93  detLayer = posSide ? tracker->posPixelForwardLayers()[layer - 1] : tracker->negPixelForwardLayers()[layer - 1];
94  }
95 
96  // Store the results for this hit.
97  Result result;
98  result.detLayer = detLayer;
99 
100  // Check that the track crosses this layer, and get the track trajectory at the crossing point.
101  std::pair<bool, TrajectoryStateOnSurface> layCross = detLayer->compatible(initTSOS, propagator, estimator);
102  if (layCross.first) {
103  LogDebug("GTTI")<<"crossed layer at "<<" x="<<layCross.second.globalPosition().x()<<" y="<<layCross.second.globalPosition().y()<<" z="<<layCross.second.globalPosition().z();
104 
105  // Find the sensor in this layer which is closest to the track trajectory.
106  // And get the track trajectory at that sensor.
107  const PropagationDirection along = alongMomentum;
108  propagator.setPropagationDirection(along);
109  std::vector< GeometricSearchDet::DetWithState > detWithState = detLayer->compatibleDets(initTSOS, propagator, estimator);
110  // Check that at least one sensor was compatible with the track trajectory.
111  if(detWithState.size() > 0) {
112  // Store track trajectory at this sensor.
113  result.valid = true;
114  result.accurate = true;
115  result.detTSOS = detWithState.front().second;
116  LogDebug("GTTI")<<" Det in this layer compatible with TSOS: subdet="<<subDet<<" layer="<<layer;
117  LogDebug("GTTI")<<" crossed sensor at x="<<result.detTSOS.globalPosition().x()<<" y="<<result.detTSOS.globalPosition().y()<<" z="<<result.detTSOS.globalPosition().z();
118 
119  } else {
120  // Track did not cross a sensor, so store approximate result from its intercept with the layer.
121  result.valid = true;
122  result.accurate = false;
123  result.detTSOS = layCross.second;
124  LogDebug("GTTI")<<" WARNING: TSOS not compatible with any det in this layer, despite having a hit in it !";
125  }
126 
127  } else {
128  // Track trajectory did not cross layer. Pathological case.
129  result.valid = false;
130  LogDebug("GTTI")<<" WARNING: track failed to cross layer, despite having a hit in hit !";
131  }
132 
133  results.push_back(result);
134  }
135  }
136 
137  return results;
138 }
#define LogDebug(id)
virtual std::pair< bool, TrajectoryStateOnSurface > compatible(const TrajectoryStateOnSurface &ts, const Propagator &, const MeasurementEstimator &) const =0
int i
Definition: DBlmapReader.cc:9
static uint32_t getLayer(uint32_t pattern)
Definition: HitPattern.h:513
static bool trackerHitFilter(uint32_t pattern)
Definition: HitPattern.h:494
T y() const
Definition: PV3DBase.h:63
GlobalPoint globalPosition() const
PropagationDirection
virtual std::vector< DetWithState > compatibleDets(const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est) const
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:139
susybsm::HSCParticleRefProd hp
Definition: classes.h:27
const DetLayer * detLayer
T z() const
Definition: PV3DBase.h:64
tuple result
Definition: query.py:137
int numberOfHits() const
Definition: HitPattern.cc:213
static uint32_t getSubStructure(uint32_t pattern)
Definition: HitPattern.h:507
edm::ESHandle< TransientTrackBuilder > trkTool_
TrajectoryStateOnSurface detTSOS
static bool validHitFilter(uint32_t pattern)
Definition: HitPattern.h:557
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:221
std::vector< Result > analyze(const edm::EventSetup &iSetup, const reco::Track &track)
const T & get() const
Definition: EventSetup.h:55
TrajectoryStateOnSurface impactPointState() const
uint32_t getHitPattern(int position) const
Definition: HitPattern.cc:144
T x() const
Definition: PV3DBase.h:62