CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/RecoParticleFlow/PFTracking/src/PFCheckHitPattern.cc

Go to the documentation of this file.
00001 #include "RecoParticleFlow/PFTracking/interface/PFCheckHitPattern.h"
00002 
00003 // To get Tracker Geometry
00004 #include "FWCore/Framework/interface/ESHandle.h"
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/TIBDetId.h"
00013 #include "DataFormats/SiStripDetId/interface/TOBDetId.h"
00014 #include "DataFormats/SiStripDetId/interface/TECDetId.h"
00015 #include "DataFormats/SiStripDetId/interface/TIDDetId.h"
00016 #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
00017 #include "DataFormats/SiPixelDetId/interface/PXFDetId.h"
00018 
00019 #include "Utilities/General/interface/CMSexception.h"
00020 
00021 #include <map>
00022 
00023 using namespace reco;
00024 using namespace std;
00025 
00026 // For a given subdetector & layer number, this static map stores the minimum and maximum
00027 // r (or z) values if it is barrel (or endcap) respectively.
00028 PFCheckHitPattern::RZrangeMap PFCheckHitPattern::rangeRorZ_;
00029 
00030 void PFCheckHitPattern::init(edm::ESHandle<TrackerGeometry> tkerGeomHandle_) {
00031 
00032   //
00033   // Note min/max radius (z) of each barrel layer (endcap disk).
00034   //
00035 
00036   geomInitDone_ = true;
00037 
00038   // Get Tracker geometry
00039   const TrackingGeometry::DetContainer& dets = tkerGeomHandle_->dets();
00040 
00041   // Loop over all modules in the Tracker.
00042   for (unsigned int i = 0; i < dets.size(); i++) {    
00043 
00044     // Get subdet and layer of this module
00045     DetInfo detInfo = this->interpretDetId(dets[i]->geographicalId());
00046     uint32_t subDet = detInfo.first;
00047 
00048     // Note r (or z) of module if barrel (or endcap).
00049     double r_or_z;
00050     if (this->barrel(subDet)) {
00051       r_or_z = dets[i]->position().perp();
00052     } else {
00053       r_or_z = fabs(dets[i]->position().z());
00054     }
00055 
00056     // Recover min/max r/z value of this layer/disk found so far.
00057     double minRZ = 999.;
00058     double maxRZ = 0.;
00059     if (rangeRorZ_.find(detInfo) != rangeRorZ_.end()) {
00060       minRZ = rangeRorZ_[detInfo].first;
00061       maxRZ = rangeRorZ_[detInfo].second;
00062     }
00063 
00064     // Update with those of this module if it exceeds them.
00065     if (minRZ > r_or_z) minRZ = r_or_z; 
00066     if (maxRZ < r_or_z) maxRZ = r_or_z;     
00067     rangeRorZ_[detInfo] = pair<double, double>(minRZ, maxRZ);
00068   }
00069 
00070 #ifdef DEBUG_CHECKHITPATTERN
00071   RZrangeMap::const_iterator d;
00072   for (d = rangeRorZ_.begin(); d != rangeRorZ_.end(); d++) {
00073     DetInfo detInfo = d->first;
00074     pair<double, double> rangeRZ = d->second;
00075   }
00076 #endif
00077 }
00078 
00079 PFCheckHitPattern::DetInfo PFCheckHitPattern::interpretDetId(DetId detId) {
00080   // Convert detId to a pair<uint32, uint32> consisting of the numbers used by HitPattern 
00081   // to identify subdetector and layer number respectively.
00082   if (detId.subdetId() == StripSubdetector::TIB) {
00083     return DetInfo( detId.subdetId() , TIBDetId(detId).layer() );
00084   } else if (detId.subdetId() == StripSubdetector::TOB) {
00085     return DetInfo( detId.subdetId() , TOBDetId(detId).layer() );
00086   } else if (detId.subdetId() == StripSubdetector::TID) {
00087     return DetInfo( detId.subdetId() , TIDDetId(detId).wheel() );
00088   } else if (detId.subdetId() == StripSubdetector::TEC) {
00089     return DetInfo( detId.subdetId() , TECDetId(detId).wheel() );
00090   } else if (detId.subdetId() == PixelSubdetector::PixelBarrel) {
00091     return DetInfo( detId.subdetId() , PXBDetId(detId).layer() );
00092   } else if (detId.subdetId() == PixelSubdetector::PixelEndcap) {
00093     return DetInfo( detId.subdetId() , PXFDetId(detId).disk() );
00094   } else {
00095     throw Genexception("Found DetId that is not in Tracker");
00096   }   
00097 }
00098 
00099 bool PFCheckHitPattern::barrel(uint32_t subDet) {
00100   // Determines if given sub-detector is in the barrel.
00101   return (subDet == StripSubdetector::TIB || subDet == StripSubdetector::TOB ||
00102           subDet == PixelSubdetector::PixelBarrel); 
00103 }
00104 
00105 
00106 pair< PFCheckHitPattern::PFTrackHitInfo, PFCheckHitPattern::PFTrackHitInfo> 
00107 PFCheckHitPattern::analyze(edm::ESHandle<TrackerGeometry> tkerGeomHandle_, 
00108                            const TrackBaseRef track, const TransientVertex& vert) 
00109 {
00110 
00111   // PFCheck if hit pattern of this track is consistent with it being produced
00112   // at given vertex. Pair.first gives number of hits on track in front of vertex.
00113   // Pair.second gives number of missing hits between vertex and innermost hit
00114   // on track.
00115 
00116   // Initialise geometry info if not yet done.
00117   if (!geomInitDone_) this->init(tkerGeomHandle_);
00118 
00119   // Get hit patterns of this track
00120   const reco::HitPattern& hp = track.get()->hitPattern(); 
00121   const reco::HitPattern& ip = track.get()->trackerExpectedHitsInner(); 
00122 
00123   // Count number of valid hits on track definately in front of the vertex,
00124   // taking into account finite depth of each layer.
00125   unsigned int nHitBefore = 0;
00126   unsigned int nHitAfter = 0;
00127 
00128   for (int i = 0; i < hp.numberOfHits(); i++) {
00129     uint32_t hit = hp.getHitPattern(i);
00130     if (hp.trackerHitFilter(hit) && hp.validHitFilter(hit)) {
00131       uint32_t subDet = hp.getSubStructure(hit);
00132       uint32_t layer = hp.getLayer(hit);
00133       DetInfo detInfo(subDet, layer);
00134       double maxRZ = rangeRorZ_[detInfo].second;
00135 
00136       if (this->barrel(subDet)) {
00137         if (vert.position().perp() > maxRZ) nHitBefore++;
00138         else nHitAfter++;
00139       } else {
00140         if (fabs(vert.position().z()) > maxRZ) nHitBefore++;
00141         else nHitAfter++;
00142       } 
00143     }
00144   }
00145 
00146   // Count number of missing hits before the innermost hit on the track,
00147   // taking into account finite depth of each layer.
00148   unsigned int nMissHitAfter = 0;
00149   unsigned int nMissHitBefore = 0;
00150 
00151   for (int i = 0; i < ip.numberOfHits(); i++) {
00152     uint32_t hit = ip.getHitPattern(i);
00153     if (ip.trackerHitFilter(hit)) {
00154       uint32_t subDet = ip.getSubStructure(hit);
00155       uint32_t layer = ip.getLayer(hit);
00156       DetInfo detInfo(subDet, layer);
00157       double minRZ = rangeRorZ_[detInfo].first;
00158 
00159       //      cout << "subDet = " << subDet << " layer = " << layer << " minRZ = " << minRZ << endl;
00160 
00161       if (this->barrel(subDet)) {
00162         if (vert.position().perp() < minRZ) nMissHitAfter++;
00163         else nMissHitBefore++;
00164       } else {
00165         if (fabs(vert.position().z()) < minRZ) nMissHitAfter++; 
00166         else nMissHitBefore++;
00167       } 
00168     }
00169   }
00170 
00171 
00172   PFTrackHitInfo trackToVertex(nHitBefore, nMissHitBefore);
00173   PFTrackHitInfo trackFromVertex(nHitAfter, nMissHitAfter);
00174 
00175 
00176   return pair< PFTrackHitInfo, PFTrackHitInfo>(trackToVertex, trackFromVertex);
00177 }
00178 
00179 void PFCheckHitPattern::print(const TrackBaseRef track) const {
00180   // Get hit patterns of this track
00181   const reco::HitPattern& hp = track.get()->hitPattern(); 
00182   const reco::HitPattern& ip = track.get()->trackerExpectedHitsInner(); 
00183 
00184   cout<<"=== Hits on Track ==="<<endl;
00185   this->print(hp);
00186   cout<<"=== Hits before track ==="<<endl;
00187   this->print(ip);
00188 }
00189 
00190 void PFCheckHitPattern::print(const reco::HitPattern& hp) const {
00191   for (int i = 0; i < hp.numberOfHits(); i++) {
00192     uint32_t hit = hp.getHitPattern(i);
00193     if (hp.trackerHitFilter(hit)) {
00194       uint32_t subdet = hp.getSubStructure(hit);
00195       uint32_t layer = hp.getLayer(hit);
00196       cout<<"hit "<<i<<" subdet="<<subdet<<" layer="<<layer<<" type "<<hp.getHitType(hit)<<endl;
00197     }
00198   } 
00199 }