CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/RecoTracker/DebugTools/src/GetTrackTrajInfo.cc

Go to the documentation of this file.
00001 #include "RecoTracker/DebugTools/interface/GetTrackTrajInfo.h"
00002 
00003 #include "RecoTracker/Record/interface/TrackerRecoGeometryRecord.h"
00004 #include "RecoTracker/TkDetLayers/interface/GeometricSearchTracker.h"
00005 #include "RecoTracker/Record/interface/CkfComponentsRecord.h"
00006 
00007 // To convert detId to subdet/layer number.
00008 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
00009 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00010 //#include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
00011 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
00012 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00013 
00014 #include "TrackingTools/KalmanUpdators/interface/Chi2MeasurementEstimator.h"
00015 #include "TrackingTools/GeomPropagators/interface/AnalyticalPropagator.h"
00016 
00017 #include "TrackingTools/DetLayers/interface/NavigationSchool.h"
00018 #include "RecoTracker/Record/interface/NavigationSchoolRecord.h"
00019 #include "TrackingTools/DetLayers/interface/NavigationSetter.h"
00020 
00021 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00022 #include "DataFormats/Math/interface/deltaPhi.h"
00023 
00024 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00025 
00026 #include "FWCore/Utilities/interface/Exception.h"
00027 
00028 std::vector< GetTrackTrajInfo::Result > GetTrackTrajInfo::analyze(const edm::EventSetup& iSetup, const reco::Track& track) 
00029 {
00030   // Determine the track trajectory and detLayer at each layer that the track produces a hit in.
00031 
00032   std::vector< GetTrackTrajInfo::Result > results;
00033 
00034   // Initialise Tracker geometry info (not sufficient to do this only on first call).
00035   edm::ESHandle<GeometricSearchTracker> tracker;
00036   iSetup.get<TrackerRecoGeometryRecord>().get( tracker );    
00037 
00038   // This is also needed to extrapolate amongst the tracker layers.
00039   edm::ESHandle<NavigationSchool> theSchool;
00040   iSetup.get<NavigationSchoolRecord>().get("SimpleNavigationSchool",theSchool);
00041   NavigationSetter junk(*theSchool);
00042 
00043   // Get the magnetic field and use it to define a propagator for extrapolating the track trajectory.
00044   edm::ESHandle<MagneticField> magField;
00045   iSetup.get<IdealMagneticFieldRecord>().get(magField);
00046   AnalyticalPropagator  propagator(&(*magField), alongMomentum);
00047 
00048   // This is used to check if a track is compatible with crossing a sensor.
00049   // Use +3.0 rather than default -3.0 here, so hit defined as inside acceptance if 
00050   // no more than 3*sigma outside detector edge, as opposed to more than 3*sigma inside detector edge.
00051   Chi2MeasurementEstimator estimator(30.,3.0);
00052 
00053   // Convert track to transientTrack, and hence to TSOS at its point of closest approach to beam.
00054   iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",trkTool_); // Needed for vertex fits
00055   reco::TransientTrack t_trk = trkTool_->build(track);
00056   TrajectoryStateOnSurface initTSOS = t_trk.impactPointState();
00057   LogDebug("GTTI")<<"TRACK TSOS POS: x="<<initTSOS.globalPosition().x()<<" y="<<initTSOS.globalPosition().y()<<" z="<<initTSOS.globalPosition().z();
00058 
00059   // Note if the track is going into +ve or -ve z.
00060   // This is only used to guess if the track is more likely to have hit a +ve rather than a -ve endcap
00061   // disk. Since the +ve and -ve disks are a long way apart, this approximate method is good enough.
00062   // More precise would be to check both possiblities and see which one (if either) the track crosses
00063   // using detLayer::compatible().
00064   bool posSide = track.eta() > 0;
00065 
00066   // Get hit patterns of this track
00067   const reco::HitPattern& hp = track.hitPattern(); 
00068 
00069   // Loop over info for each hit
00070   // N.B. Hits are sorted according to increasing distance from origin by
00071   // RecoTracker/TrackProducer/src/TrackProducerBase.cc
00072   for (int i = 0; i < hp.numberOfHits(); i++) {
00073     uint32_t hit = hp.getHitPattern(i);
00074     if (hp.trackerHitFilter(hit) && hp.validHitFilter(hit)) {
00075       uint32_t subDet = hp.getSubStructure(hit);
00076       uint32_t layer = hp.getLayer(hit);
00077       // subdet: PixelBarrel=1, PixelEndcap=2, TIB=3, TID=4, TOB=5, TEC=6
00078       LogDebug("GTTI")<<"    hit in subdet="<<subDet<<" layer="<<layer;
00079 
00080       // Get corresponding DetLayer object (based on code in GeometricSearchTracker::idToLayer(...)
00081       const DetLayer* detLayer = 0;
00082       if (subDet == StripSubdetector::TIB) {
00083         detLayer = tracker->tibLayers()[layer - 1];
00084       } else if (subDet == StripSubdetector::TOB) {
00085         detLayer = tracker->tobLayers()[layer - 1];
00086       } else if (subDet == StripSubdetector::TID) {
00087         detLayer = posSide ? tracker->posTidLayers()[layer - 1] : tracker->negTidLayers()[layer - 1];
00088       } else if (subDet == StripSubdetector::TEC) {
00089         detLayer = posSide ? tracker->posTecLayers()[layer - 1] : tracker->negTecLayers()[layer - 1];
00090       } else if (subDet == PixelSubdetector::PixelBarrel) {
00091         detLayer = tracker->pixelBarrelLayers()[layer - 1];
00092       } else if (subDet == PixelSubdetector::PixelEndcap) {
00093         detLayer = posSide ? tracker->posPixelForwardLayers()[layer - 1] : tracker->negPixelForwardLayers()[layer - 1];
00094       }
00095 
00096       // Store the results for this hit.
00097       Result result;
00098       result.detLayer = detLayer;
00099 
00100       // Check that the track crosses this layer, and get the track trajectory at the crossing point.
00101       std::pair<bool, TrajectoryStateOnSurface> layCross = detLayer->compatible(initTSOS, propagator, estimator);
00102       if (layCross.first) {
00103         LogDebug("GTTI")<<"crossed layer at "<<" x="<<layCross.second.globalPosition().x()<<" y="<<layCross.second.globalPosition().y()<<" z="<<layCross.second.globalPosition().z();
00104 
00105         // Find the sensor in this layer which is closest to the track trajectory.
00106         // And get the track trajectory at that sensor.
00107         const PropagationDirection along = alongMomentum;
00108         propagator.setPropagationDirection(along);
00109         std::vector< GeometricSearchDet::DetWithState > detWithState = detLayer->compatibleDets(initTSOS, propagator, estimator);
00110         // Check that at least one sensor was compatible with the track trajectory.
00111         if(detWithState.size() > 0) {
00112           // Store track trajectory at this sensor.
00113           result.valid    = true;
00114           result.accurate = true;
00115           result.detTSOS  = detWithState.front().second;
00116           LogDebug("GTTI")<<"      Det in this layer compatible with TSOS: subdet="<<subDet<<" layer="<<layer;
00117           LogDebug("GTTI")<<"      crossed sensor at x="<<result.detTSOS.globalPosition().x()<<" y="<<result.detTSOS.globalPosition().y()<<" z="<<result.detTSOS.globalPosition().z();
00118 
00119         } else {
00120           // Track did not cross a sensor, so store approximate result from its intercept with the layer.
00121           result.valid    = true;
00122           result.accurate = false;
00123           result.detTSOS  = layCross.second;
00124           LogDebug("GTTI")<<"      WARNING: TSOS not compatible with any det in this layer, despite having a hit in it !";
00125         }
00126 
00127       } else {
00128         // Track trajectory did not cross layer. Pathological case.
00129         result.valid = false;
00130         LogDebug("GTTI")<<"      WARNING: track failed to cross layer, despite having a hit in hit !";
00131       }
00132 
00133       results.push_back(result);
00134     }
00135   }
00136 
00137   return results;
00138 }