CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoPixelVertexing/PixelLowPtUtilities/src/ClusterShapeTrackFilter.cc

Go to the documentation of this file.
00001 #include "RecoPixelVertexing/PixelLowPtUtilities/interface/ClusterShapeTrackFilter.h"
00002 
00003 #include "RecoPixelVertexing/PixelTrackFitting/src/CircleFromThreePoints.h"
00004 #include "RecoPixelVertexing/PixelLowPtUtilities/interface/HitInfo.h"
00005 #include "RecoPixelVertexing/PixelLowPtUtilities/interface/ClusterShapeHitFilter.h"
00006 
00007 #include "FWCore/Framework/interface/EventSetup.h"
00008 #include "FWCore/Framework/interface/ESHandle.h"
00009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00011 
00012 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
00013 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
00014 #include "DataFormats/TrackReco/interface/Track.h"
00015 
00016 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00017 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00018 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
00019 
00020 #include "RecoTracker/Record/interface/CkfComponentsRecord.h"
00021 
00022 inline float sqr(float x) { return x*x; }
00023 
00024 using namespace std;
00025 
00026 /*****************************************************************************/
00027 ClusterShapeTrackFilter::ClusterShapeTrackFilter
00028   (const edm::ParameterSet& ps, const edm::EventSetup& es)
00029 {
00030   // Get tracker geometry
00031   edm::ESHandle<TrackerGeometry> tracker;
00032   es.get<TrackerDigiGeometryRecord>().get(tracker);
00033   theTracker = tracker.product();
00034 
00035   // Get cluster shape hit filter
00036   edm::ESHandle<ClusterShapeHitFilter> shape;
00037   es.get<CkfComponentsRecord>().get("ClusterShapeHitFilter",shape);
00038   theFilter = shape.product();
00039 
00040   // Get ptMin if available
00041   ptMin = (ps.exists("ptMin") ? ps.getParameter<double>("ptMin") : 0.);
00042   ptMax = (ps.exists("ptMax") ? ps.getParameter<double>("ptMax") : 999999.);
00043 }
00044 
00045 /*****************************************************************************/
00046 ClusterShapeTrackFilter::~ClusterShapeTrackFilter()
00047 {
00048 }
00049 
00050 /*****************************************************************************/
00051 float ClusterShapeTrackFilter::areaParallelogram
00052   (const Global2DVector& a, const Global2DVector& b) const
00053 {  
00054   return a.x() * b.y() - a.y() * b.x();
00055 }
00056 
00057 /*****************************************************************************/
00058 vector<GlobalVector> ClusterShapeTrackFilter::getGlobalDirs
00059   (const vector<GlobalPoint> & g) const
00060 {
00061   // Get 2d points
00062   vector<Global2DVector> p;
00063   for(vector<GlobalPoint>::const_iterator ig = g.begin();
00064                                           ig!= g.end(); ig++)
00065      p.push_back( Global2DVector(ig->x(), ig->y()) );
00066 
00067   //
00068   vector<GlobalVector> globalDirs;
00069 
00070   // Determine circle
00071   CircleFromThreePoints circle(g[0],g[1],g[2]);
00072 
00073   if(circle.curvature() != 0.)
00074   {
00075     Global2DVector c (circle.center().x(), circle.center().y());
00076 
00077     float rad2 = (p[0] - c).mag2();
00078     float a12 = asin(fabsf(areaParallelogram(p[0] - c, p[1] - c)) / rad2);
00079 
00080     float slope = (g[1].z() - g[0].z()) / a12;
00081 
00082     float cotTheta = slope * circle.curvature(); // == sinhEta
00083     float coshEta  = sqrt(1 + sqr(cotTheta));    // == 1/sinTheta
00084 
00085     // Calculate globalDirs
00086     float sinTheta =       1. / coshEta;
00087     float cosTheta = cotTheta * sinTheta;
00088 
00089     int dir;
00090     if(areaParallelogram(p[0] - c, p[1] - c) > 0) dir = 1; else dir = -1;
00091 
00092     float curvature = circle.curvature();
00093 
00094     for(vector<Global2DVector>::const_iterator ip = p.begin();
00095                                                ip!= p.end(); ip++)
00096     {
00097       Global2DVector v = (*ip - c)*curvature*dir;
00098       globalDirs.push_back(GlobalVector(-v.y()*sinTheta,
00099                                          v.x()*sinTheta,
00100                                                cosTheta));
00101     }
00102   }
00103 
00104   return globalDirs;
00105 }
00106 
00107 /*****************************************************************************/
00108 vector<GlobalPoint> ClusterShapeTrackFilter::getGlobalPoss
00109   (const vector<const TrackingRecHit *> & recHits) const
00110 {
00111   vector<GlobalPoint> globalPoss;
00112 
00113   for(vector<const TrackingRecHit *>::const_iterator recHit = recHits.begin();
00114                                                      recHit!= recHits.end();
00115                                                      recHit++)
00116   {
00117     DetId detId = (*recHit)->geographicalId();
00118 
00119     GlobalPoint gpos = 
00120       theTracker->idToDet(detId)->toGlobal((*recHit)->localPosition());
00121 
00122     globalPoss.push_back(gpos);
00123   }
00124 
00125   return globalPoss;
00126 }
00127 
00128 /*****************************************************************************/
00129 bool ClusterShapeTrackFilter::operator()
00130   (const reco::Track* track,
00131    const vector<const TrackingRecHit *> & recHits) const
00132 {
00133   // Do not even look at pairs
00134   if(recHits.size() <= 2) return true;
00135 
00136   // Check pt
00137   if(track->pt() < ptMin ||
00138      track->pt() > ptMax)
00139   {
00140     LogTrace("ClusterShapeTrackFilter")
00141        << "  [ClusterShapeTrackFilter] pt not in range: "
00142        << ptMin << " " << track->pt() << " " << ptMax;
00143     return false;
00144   }
00145 
00146   // Get global positions
00147   vector<GlobalPoint>  globalPoss = getGlobalPoss(recHits);
00148 
00149   // Get global directions
00150   vector<GlobalVector> globalDirs = getGlobalDirs(globalPoss);
00151 
00152   bool ok = true;
00153 
00154   // Check whether shape of pixel cluster is compatible
00155   // with local track direction
00156   for(unsigned int i = 0; i < recHits.size(); i++)
00157   {
00158     const SiPixelRecHit* pixelRecHit =
00159       dynamic_cast<const SiPixelRecHit *>(recHits[i]);
00160 
00161     if(!pixelRecHit->isValid())
00162     { 
00163       ok = false; break; 
00164     }
00165 
00166     if(! theFilter->isCompatible(*pixelRecHit, globalDirs[i]) )
00167     {
00168       LogTrace("ClusterShapeTrackFilter")
00169          << "  [ClusterShapeTrackFilter] clusShape problem"
00170          << HitInfo::getInfo(*recHits[i]);
00171 
00172       ok = false; break;
00173     }
00174   }
00175 
00176   return ok;
00177 }
00178