CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/PhysicsTools/RecoUtils/src/CheckHitPattern.cc

Go to the documentation of this file.
00001 #include "PhysicsTools/RecoUtils/interface/CheckHitPattern.h"
00002 #include "RecoTracker/DebugTools/interface/FixTrackHitPattern.h"
00003 
00004 // To get Tracker Geometry
00005 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00006 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00007 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
00008 
00009 // To convert detId to subdet/layer number.
00010 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
00011 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00012 //#include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
00013 #include "DataFormats/SiStripDetId/interface/TIBDetId.h"
00014 #include "DataFormats/SiStripDetId/interface/TOBDetId.h"
00015 #include "DataFormats/SiStripDetId/interface/TECDetId.h"
00016 #include "DataFormats/SiStripDetId/interface/TIDDetId.h"
00017 #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
00018 #include "DataFormats/SiPixelDetId/interface/PXFDetId.h"
00019 
00020 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00021 #include "DataFormats/Math/interface/deltaPhi.h"
00022 
00023 #include "FWCore/Utilities/interface/Exception.h"
00024 
00025 // For a given subdetector & layer number, this static map stores the minimum and maximum
00026 // r (or z) values if it is barrel (or endcap) respectively.
00027 CheckHitPattern::RZrangeMap CheckHitPattern::rangeRorZ_;
00028 
00029 void CheckHitPattern::init(const edm::EventSetup& iSetup) {
00030 
00031   //
00032   // Note min/max radius (z) of each barrel layer (endcap disk).
00033   //
00034 
00035   geomInitDone_ = true;
00036 
00037   // Get Tracker geometry
00038   edm::ESHandle<TrackerGeometry> trackerGeometry;
00039   iSetup.get<TrackerDigiGeometryRecord>().get(trackerGeometry);
00040   const TrackingGeometry::DetContainer& dets = trackerGeometry->dets();
00041 
00042   // Loop over all modules in the Tracker.
00043   for (unsigned int i = 0; i < dets.size(); i++) {    
00044 
00045     // Get subdet and layer of this module
00046     DetInfo detInfo = this->interpretDetId(dets[i]->geographicalId());
00047     uint32_t subDet = detInfo.first;
00048 
00049     // Note r (or z) of module if barrel (or endcap).
00050     double r_or_z;
00051     if (this->barrel(subDet)) {
00052       r_or_z = dets[i]->position().perp();
00053     } else {
00054       r_or_z = fabs(dets[i]->position().z());
00055     }
00056 
00057     // Recover min/max r/z value of this layer/disk found so far.
00058     double minRZ = 999.;
00059     double maxRZ = 0.;
00060     if (rangeRorZ_.find(detInfo) != rangeRorZ_.end()) {
00061       minRZ = rangeRorZ_[detInfo].first;
00062       maxRZ = rangeRorZ_[detInfo].second;
00063     }
00064 
00065     // Update with those of this module if it exceeds them.
00066     if (minRZ > r_or_z) minRZ = r_or_z; 
00067     if (maxRZ < r_or_z) maxRZ = r_or_z;     
00068     rangeRorZ_[detInfo] = std::pair<double, double>(minRZ, maxRZ);
00069   }
00070 
00071 #ifdef DEBUG_CHECKHITPATTERN
00072   RZrangeMap::const_iterator d;
00073   for (d = rangeRorZ_.begin(); d != rangeRorZ_.end(); d++) {
00074     DetInfo detInfo = d->first;
00075     std::pair<double, double> rangeRZ = d->second;
00076     std::std::cout<<"CHECKHITPATTERN: Tracker subdetector type="<<detInfo.first<<" layer="<<detInfo.second
00077         <<" has min r (or z) ="<<rangeRZ.first<<" and max r (or z) = "<<rangeRZ.second<<std::std::endl; 
00078   }
00079 #endif
00080 }
00081 
00082 CheckHitPattern::DetInfo CheckHitPattern::interpretDetId(DetId detId) {
00083   // Convert detId to a pair<uint32, uint32> consisting of the numbers used by HitPattern 
00084   // to identify subdetector and layer number respectively.
00085   if (detId.subdetId() == StripSubdetector::TIB) {
00086     return DetInfo( detId.subdetId() , TIBDetId(detId).layer() );
00087   } else if (detId.subdetId() == StripSubdetector::TOB) {
00088     return DetInfo( detId.subdetId() , TOBDetId(detId).layer() );
00089   } else if (detId.subdetId() == StripSubdetector::TID) {
00090     return DetInfo( detId.subdetId() , TIDDetId(detId).wheel() );
00091   } else if (detId.subdetId() == StripSubdetector::TEC) {
00092     return DetInfo( detId.subdetId() , TECDetId(detId).wheel() );
00093   } else if (detId.subdetId() == PixelSubdetector::PixelBarrel) {
00094     return DetInfo( detId.subdetId() , PXBDetId(detId).layer() );
00095   } else if (detId.subdetId() == PixelSubdetector::PixelEndcap) {
00096     return DetInfo( detId.subdetId() , PXFDetId(detId).disk() );
00097   } else {
00098     throw cms::Exception("NotFound","Found DetId that is not in Tracker");
00099   }   
00100 }
00101 
00102 bool CheckHitPattern::barrel(uint32_t subDet) {
00103   // Determines if given sub-detector is in the barrel.
00104   return (subDet == StripSubdetector::TIB || subDet == StripSubdetector::TOB ||
00105           subDet == PixelSubdetector::PixelBarrel); 
00106 }
00107 
00108 
00109 CheckHitPattern::Result CheckHitPattern::analyze(const edm::EventSetup& iSetup, 
00110                          const reco::Track& track, const VertexState& vert, bool fixHitPattern) 
00111 {
00112   // Check if hit pattern of this track is consistent with it being produced
00113   // at given vertex. 
00114 
00115   // Initialise geometry info if not yet done.
00116   if (!geomInitDone_) this->init(iSetup);
00117 
00118   // Optionally set vertex position to zero for debugging.
00119   // VertexState vertDebug( GlobalPoint(0.,0.,0.) , GlobalError(1e-8, 0., 1e-8, 0., 0., 1e-8) );
00120 
00121   // Evaluate track parameters at vertex.
00122   iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",trkTool_); // Needed for vertex fits
00123   reco::TransientTrack t_trk = trkTool_->build(track);
00124   GlobalVector p3_trk = t_trk.trajectoryStateClosestToPoint(vert.position()).momentum();
00125   bool trkGoesInsideOut = fabs(reco::deltaPhi<const GlobalVector, const GlobalPoint>(p3_trk, vert.position())) < 0.5*M_PI;
00126 
00127   LogDebug("CHP")<<"TRACK: in-->out ? "<<trkGoesInsideOut<<" dxy="<<track.dxy()<<" sz="<<track.dz()<<" eta="<<track.eta()<<" barrel hits="<<track.hitPattern().numberOfValidPixelHits()<<"/"<<track.hitPattern().numberOfValidStripTIBHits()<<"/"<<track.hitPattern().numberOfValidStripTOBHits();
00128   LogDebug("CHP")<<"VERT: r="<<vert.position().perp()<<" z="<<vert.position().z();
00129   //  if (vert.position().perp() < 3.5 && fabs(vert.position().z()) < 10. && fabs(track.eta()) < 1 && fabs(track.dxy()) < 2 && fabs(track.dz()) < 2 && track.hitPattern().numberOfValidPixelHits() == 0 && track.hitPattern().numberOfValidStripTIBHits() == 0) LogDebug("CHP")<<"LOOKATTHISTRACK";
00130   // Get hit patterns of this track
00131   const reco::HitPattern& hp = track.hitPattern(); 
00132   reco::HitPattern        ip = track.trackerExpectedHitsInner(); 
00133 
00134   // Optionally fix inner hit pattern (needed if uncertainty on track trajectory is large).
00135   if (fixHitPattern) {
00136     static FixTrackHitPattern fixTrackHitPattern;
00137     FixTrackHitPattern::Result fixedHP = fixTrackHitPattern.analyze(iSetup, track);
00138     ip = fixedHP.innerHitPattern;
00139   }
00140   
00141   // Count number of valid hits on track definately in front of the vertex,
00142   // taking into account finite depth of each layer.
00143   unsigned int nHitBefore = 0;
00144   for (int i = 0; i < hp.numberOfHits(); i++) {
00145     uint32_t hit = hp.getHitPattern(i);
00146     if (hp.trackerHitFilter(hit) && hp.validHitFilter(hit)) {
00147       uint32_t subDet = hp.getSubStructure(hit);
00148       uint32_t layer = hp.getLayer(hit);
00149       DetInfo detInfo(subDet, layer);
00150       double maxRZ = rangeRorZ_[detInfo].second;
00151 
00152       if (this->barrel(subDet)) {
00153         // Be careful. If the track starts by going outside-->in, it is allowed to have hits before the vertex !
00154         if (vert.position().perp() > maxRZ && trkGoesInsideOut) nHitBefore++;
00155       } else {
00156         if (fabs(vert.position().z()) > maxRZ) nHitBefore++;
00157       } 
00158     }
00159   }
00160 
00161   // Count number of missing hits before the innermost hit on the track,
00162   // taking into account finite depth of each layer.
00163   unsigned int nMissHitAfter = 0;
00164   for (int i = 0; i < ip.numberOfHits(); i++) {
00165     uint32_t hit = ip.getHitPattern(i);
00166     //    if (ip.trackerHitFilter(hit)) {
00167     if (ip.trackerHitFilter(hit) && ip.type_1_HitFilter(hit)) {
00168       uint32_t subDet = ip.getSubStructure(hit);
00169       uint32_t layer = ip.getLayer(hit);
00170       DetInfo detInfo(subDet, layer);
00171       double minRZ = rangeRorZ_[detInfo].first;
00172 
00173       if (this->barrel(subDet)) {
00174         // Be careful. If the track starts by going outside-->in, then it misses hits
00175         // in all layers it crosses  before its innermost valid hit.
00176         if (vert.position().perp() < minRZ || ! trkGoesInsideOut) nMissHitAfter++;
00177       } else {
00178         if (fabs(vert.position().z()) < minRZ) nMissHitAfter++;
00179       } 
00180     }
00181   }
00182  
00183   Result result;
00184   result.hitsInFrontOfVert = nHitBefore;
00185   result.missHitsAfterVert = nMissHitAfter;
00186   return result;
00187 }
00188 
00189 void CheckHitPattern::print(const reco::Track& track) const {
00190   // Get hit patterns of this track
00191   const reco::HitPattern& hp = track.hitPattern(); 
00192   const reco::HitPattern& ip = track.trackerExpectedHitsInner(); 
00193 
00194   std::cout<<"=== Hits on Track ==="<<std::endl;
00195   this->print(hp);
00196   std::cout<<"=== Hits before track ==="<<std::endl;
00197   this->print(ip);
00198 }
00199 
00200 void CheckHitPattern::print(const reco::HitPattern& hp) const {
00201   for (int i = 0; i < hp.numberOfHits(); i++) {
00202     uint32_t hit = hp.getHitPattern(i);
00203     if (hp.trackerHitFilter(hit)) {
00204       uint32_t subdet = hp.getSubStructure(hit);
00205       uint32_t layer = hp.getLayer(hit);
00206       std::cout<<"hit "<<i<<" subdet="<<subdet<<" layer="<<layer<<" type "<<hp.getHitType(hit)<<std::endl;
00207     }
00208   } 
00209 }