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