CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/RecoTracker/DebugTools/src/FixTrackHitPattern.cc

Go to the documentation of this file.
00001 #include "RecoTracker/DebugTools/interface/FixTrackHitPattern.h"
00002 #include "RecoTracker/DebugTools/interface/GetTrackTrajInfo.h"
00003 
00004 #include "RecoTracker/Record/interface/TrackerRecoGeometryRecord.h"
00005 #include "RecoTracker/TkDetLayers/interface/GeometricSearchTracker.h"
00006 #include "RecoTracker/MeasurementDet/interface/MeasurementTracker.h"
00007 #include "TrackingTools/MeasurementDet/interface/MeasurementDet.h"
00008 #include "RecoTracker/Record/interface/CkfComponentsRecord.h"
00009 #include "DataFormats/TrackingRecHit/interface/InvalidTrackingRecHit.h"
00010 
00011 // To convert detId to subdet/layer number.
00012 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
00013 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00014 //#include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
00015 #include "DataFormats/SiStripDetId/interface/TIBDetId.h"
00016 #include "DataFormats/SiStripDetId/interface/TOBDetId.h"
00017 #include "DataFormats/SiStripDetId/interface/TECDetId.h"
00018 #include "DataFormats/SiStripDetId/interface/TIDDetId.h"
00019 #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
00020 #include "DataFormats/SiPixelDetId/interface/PXFDetId.h"
00021 
00022 #include "TrackingTools/KalmanUpdators/interface/Chi2MeasurementEstimator.h"
00023 #include "TrackingTools/GeomPropagators/interface/AnalyticalPropagator.h"
00024 
00025 #include "TrackingTools/DetLayers/interface/NavigationSchool.h"
00026 #include "RecoTracker/Record/interface/NavigationSchoolRecord.h"
00027 #include "TrackingTools/DetLayers/interface/NavigationSetter.h"
00028 
00029 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00030 #include "DataFormats/Math/interface/deltaPhi.h"
00031 
00032 #include "FWCore/Utilities/interface/Exception.h"
00033 
00034 FixTrackHitPattern::Result FixTrackHitPattern::analyze(const edm::EventSetup& iSetup, const reco::Track& track) 
00035 {
00036   // Recalculate the inner and outer missing hit patterns. See header file for detailed comments.
00037 
00038   Result result;
00039 
00040   using namespace std;
00041 
00042   // Initialise Tracker geometry info (not sufficient to do this only on first call).
00043   edm::ESHandle<GeometricSearchTracker> tracker;
00044   iSetup.get<TrackerRecoGeometryRecord>().get( tracker );    
00045 
00046   // This is also needed to extrapolate amongst the tracker layers.
00047   edm::ESHandle<NavigationSchool> theSchool;
00048   iSetup.get<NavigationSchoolRecord>().get("SimpleNavigationSchool",theSchool);
00049   NavigationSetter junk(*theSchool);
00050 
00051   // This is needed to determine which sensors are functioning.
00052   edm::ESHandle<MeasurementTracker> measTk;
00053   iSetup.get<CkfComponentsRecord>().get(measTk);
00054   // Don't do this. It tries getting the tracker clusters, which do not exist in AOD.
00055   // Hopefully not needed if one simply wants to know which sensors are active.
00056   // measTk->update(iEvent);
00057 
00058   // Get the magnetic field and use it to define a propagator for extrapolating the track trajectory.
00059   edm::ESHandle<MagneticField> magField;
00060   iSetup.get<IdealMagneticFieldRecord>().get(magField);
00061   AnalyticalPropagator  propagator(&(*magField), alongMomentum);
00062 
00063   // This is used to check if a track is compatible with crossing a sensor.
00064   // Use +3.0 rather than default -3.0 here, so hit defined as inside acceptance if 
00065   // no more than 3*sigma outside detector edge, as opposed to more than 3*sigma inside detector edge.
00066   Chi2MeasurementEstimator estimator(30.,3.0);
00067 
00068   // Get the track trajectory and detLayer at each valid hit on the track.
00069   static GetTrackTrajInfo getTrackTrajInfo;
00070   vector<GetTrackTrajInfo::Result> trkTrajInfo = getTrackTrajInfo.analyze(iSetup, track);
00071   unsigned int nValidHits = trkTrajInfo.size();
00072 
00073   // Calculate the inner and outer hitPatterns on the track.
00074   // i.e. The list of where the track should have had hits inside its innermost valid hit 
00075   // and outside its outmost valid hit.
00076 
00077   enum InnerOuter {INNER = 1, OUTER=2};
00078   for (unsigned int inOut = INNER; inOut <= OUTER; inOut++) {
00079 
00080     // Get the track trajectory and detLayer at the inner/outermost valid hit.
00081     unsigned int ihit = (inOut == INNER) ? 0 : nValidHits - 1;
00082 
00083     // Check that information about the track trajectory was available for this hit.
00084     if (trkTrajInfo[ihit].valid) {
00085       const DetLayer* detLayer = trkTrajInfo[ihit].detLayer;
00086       const TrajectoryStateOnSurface& detTSOS = trkTrajInfo[ihit].detTSOS;
00087       const FreeTrajectoryState* detFTS = detTSOS.freeTrajectoryState();
00088 
00089       // When getting inner hit pattern, must propagate track inwards from innermost layer with valid hit.
00090       // When getting outer hit pattern, must propagate track outwards from outermost layer with valid hit.
00091       const PropagationDirection direc = (inOut == INNER) ? oppositeToMomentum : alongMomentum;
00092 
00093       // Find all layers this track is compatible with in the desired direction, starting from this layer.
00094       // Based on code in RecoTracker/TrackProducer/interface/TrackProducerBase.icc
00095       // N.B. The following call uses code in RecoTracker/TkNavigation/src/SimpleBarrelNavigableLayer::compatibleLayers() 
00096       // and RecoTracker/TkNavigation/src/SimpleNavigableLayer::wellInside(). 
00097       // These make some curious checks on the direction of the trajectory relative to its starting point,
00098       // so care was required above when calculating detFTS.
00099       vector<const DetLayer*> compLayers = detLayer->compatibleLayers(*detFTS, direc);
00100       LogDebug("FTHP")<<"Number of inner/outer "<<inOut<<" layers intercepted by track = "<<compLayers.size()<<endl;
00101 
00102       int counter = 0;
00103       reco::HitPattern newHitPattern;
00104 
00105       for(vector<const DetLayer *>::const_iterator it=compLayers.begin(); it!=compLayers.end();
00106           ++it){
00107         if ((*it)->basicComponents().empty()) {
00108           //this should never happen. but better protect for it
00109           edm::LogWarning("FixTrackHitPattern")<<"a detlayer with no components: I can not figure out a DetId from this layer. please investigate.";
00110           continue;
00111         }
00112 
00113         // Find the sensor in this layer which is closest to the track trajectory.
00114         // And get the track trajectory at that sensor.
00115         propagator.setPropagationDirection(direc);
00116         vector< GeometricSearchDet::DetWithState > detWithState = (*it)->compatibleDets(detTSOS, propagator, estimator);
00117         // Check that at least one sensor was compatible with the track trajectory.
00118         if(detWithState.size() > 0) {
00119           // Check that this sensor is functional
00120           DetId id = detWithState.front().first->geographicalId();
00121           const MeasurementDet* measDet = measTk->idToDet(id);      
00122           if(measDet->isActive()){    
00123             // Hence record that the track should have produced a hit here, but did not.
00124             // Store the information in a HitPattern.
00125             InvalidTrackingRecHit  tmpHit(id, TrackingRecHit::missing);
00126             newHitPattern.set(tmpHit, counter);      
00127             counter++; 
00128           } else {
00129             // Missing hit expected here, since sensor was not functioning.
00130           }
00131         }
00132       }//end loop over compatible layers
00133 
00134       // Store this result.
00135       if (inOut == INNER) {
00136         result.innerHitPattern = newHitPattern;
00137       } else {
00138         result.outerHitPattern = newHitPattern;
00139       }
00140 
00141       // Print result for debugging.
00142       LogDebug("FTHP")<<"Number of missing hits "<<newHitPattern.numberOfHits()<<"/"<<counter<<endl;
00143       for (int j = 0; j < std::max(newHitPattern.numberOfHits(), counter); j++) {
00144         uint32_t hp = newHitPattern.getHitPattern(j);
00145         uint32_t subDet = newHitPattern.getSubStructure(hp);
00146         uint32_t layer = newHitPattern.getLayer(hp);
00147         uint32_t status = newHitPattern.getHitType(hp);
00148         LogDebug("FTHP")<<"           layer with no matched hit at counter="<<j<<" subdet="<<subDet<<" layer="<<layer<<" status="<<status<<endl;
00149       }
00150 
00151     } else {
00152       LogDebug("FTHP")<<"WARNING: could not calculate inner/outer hit pattern as trajectory info for inner/out hit missing"<<endl;
00153     }
00154 
00155   }
00156 
00157   return result;
00158 }