CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch2/src/RecoPixelVertexing/PixelLowPtUtilities/src/LowPtClusterShapeSeedComparitor.cc

Go to the documentation of this file.
00001 #include "RecoPixelVertexing/PixelLowPtUtilities/interface/LowPtClusterShapeSeedComparitor.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 "RecoTracker/TkSeedingLayers/interface/SeedingHitSet.h"
00015 
00016 #include "RecoTracker/Record/interface/CkfComponentsRecord.h"
00017 
00018 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00019 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
00020 #include "DataFormats/GeometryVector/interface/Basic2DVector.h"
00021 
00022 #include<cmath>
00023 
00024 
00025 
00026 namespace {
00027   typedef Basic2DVector<float>   Vector2D;
00028 
00029   inline float sqr(float x) { return x*x; }
00030 
00031   /*****************************************************************************/
00032   inline
00033   float areaParallelogram
00034   (const Vector2D& a, const Vector2D& b)
00035   {  
00036     return a.x() * b.y() - a.y() * b.x();
00037   }
00038   
00039   /*****************************************************************************/
00040 
00041   inline 
00042   bool getGlobalDirs(GlobalPoint const * g, GlobalVector * globalDirs)
00043   {
00044     
00045     
00046     // Determine circle
00047     CircleFromThreePoints circle(g[0],g[1],g[2]);
00048     
00049     float curvature = circle.curvature();
00050     if(0.f == curvature) {
00051       LogDebug("LowPtClusterShapeSeedComparitor")<<"the curvature is null:"
00052                                                  <<"\n point1: "<<g[0]
00053                                                  <<"\n point2: "<<g[1]
00054                                                  <<"\n point3: "<<g[2];
00055       return false;
00056     }
00057 
00058    // Get 2d points
00059     Vector2D p[3];
00060     Vector2D c  = circle.center();
00061     for(int i=0; i!=3; i++)
00062       p[i] =  g[i].basicVector().xy() -c;
00063  
00064 
00065     float area = std::abs(areaParallelogram(p[1] - p[0], p[1]));
00066     
00067     float a12 = std::asin(std::min(area*curvature*curvature,1.f));
00068     
00069     float slope = (g[1].z() - g[0].z()) / a12;
00070  
00071     // Calculate globalDirs
00072    
00073     float cotTheta = slope * curvature; 
00074     float sinTheta = 1.f/std::sqrt(1.f + sqr(cotTheta));
00075     float cosTheta = cotTheta*sinTheta;
00076     
00077     if (areaParallelogram(p[0], p[1] ) < 0)  sinTheta = - sinTheta;
00078         
00079     for(int i = 0; i!=3;  i++) {
00080       Vector2D vl = p[i]*(curvature*sinTheta);
00081       globalDirs[i] = GlobalVector(-vl.y(),
00082                                     vl.x(),
00083                                     cosTheta
00084                                     );
00085     }
00086     return true;
00087   }
00088 
00089   /*****************************************************************************/
00090 
00091   
00092   inline
00093   void getGlobalPos(const SeedingHitSet &hits, GlobalPoint * globalPoss)
00094   {
00095     
00096     for(unsigned int i=0; i!=hits.size(); ++i)
00097         globalPoss[i] = hits[i]->globalPosition();
00098   }
00099 
00100 } // namespace
00101 
00102 /*****************************************************************************/
00103 void LowPtClusterShapeSeedComparitor::init(const edm::EventSetup& es) {
00104   es.get<CkfComponentsRecord>().get("ClusterShapeHitFilter", theShapeFilter);
00105 }
00106 
00107 bool LowPtClusterShapeSeedComparitor::compatible(const SeedingHitSet &hits, const TrackingRegion &) const
00108 //(const reco::Track* track, const vector<const TrackingRecHit *> & recHits) const
00109 {
00110   assert(hits.size()==3);
00111 
00112   const ClusterShapeHitFilter * filter = theShapeFilter.product();
00113   assert(filter != 0 && "LowPtClusterShapeSeedComparitor: init(EventSetup) method was not called");
00114 
00115    // Get global positions
00116    GlobalPoint  globalPoss[3];
00117    getGlobalPos(hits, globalPoss);
00118 
00119   // Get global directions
00120   GlobalVector globalDirs[3]; 
00121 
00122   bool ok = getGlobalDirs(globalPoss,globalDirs);
00123 
00124   // Check whether shape of pixel cluster is compatible
00125   // with local track direction
00126 
00127   if (!ok)
00128     {
00129       LogDebug("LowPtClusterShapeSeedComparitor")<<"curvarture 0:"
00130                                                  <<"\nnHits: "<<hits.size()
00131                                                  <<" will say the seed is good anyway.";
00132       return true;
00133     }
00134 
00135   for(int i = 0; i < 3; i++)
00136   {
00137     const SiPixelRecHit* pixelRecHit =
00138       dynamic_cast<const SiPixelRecHit *>(hits[i]->hit());
00139 
00140     if (!pixelRecHit){
00141       edm::LogError("LowPtClusterShapeSeedComparitor")<<"this is not a pixel cluster";
00142       ok = false; break;
00143     }
00144 
00145     if(!pixelRecHit->isValid())
00146     { 
00147       ok = false; break; 
00148     }
00149     
00150     LogDebug("LowPtClusterShapeSeedComparitor")<<"about to compute compatibility."
00151                                                <<"hit ptr: "<<pixelRecHit
00152                                                <<"global direction:"<< globalDirs[i];
00153 
00154 
00155     if(! filter->isCompatible(*pixelRecHit, globalDirs[i]) )
00156     {
00157       LogTrace("LowPtClusterShapeSeedComparitor")
00158          << " clusShape is not compatible"
00159          << HitInfo::getInfo(*hits[i]->hit());
00160 
00161       ok = false; break;
00162     }
00163   }
00164 
00165   return ok;
00166 }
00167