CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_9_patch3/src/RecoPixelVertexing/PixelLowPtUtilities/plugins/ClusterShapeSeedComparitor.cc

Go to the documentation of this file.
00001 #include "RecoTracker/TkSeedingLayers/interface/SeedComparitor.h"
00002 #include "RecoTracker/TkSeedingLayers/interface/SeedComparitorFactory.h"
00003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00004 #include "RecoPixelVertexing/PixelLowPtUtilities/interface/ClusterShapeHitFilter.h"
00005 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
00006 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
00007 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
00008 #include "DataFormats/TrackerRecHit2D/interface/ProjectedSiStripRecHit2D.h"
00009 #include "FWCore/Framework/interface/ESHandle.h"
00010 #include "RecoTracker/Record/interface/CkfComponentsRecord.h"
00011 #include "RecoTracker/TkSeedGenerator/interface/FastHelix.h"
00012 #include "RecoTracker/TkSeedingLayers/interface/SeedingHitSet.h"
00013 #include "MagneticField/Engine/interface/MagneticField.h"
00014 #include <cstdio>
00015 #include <cassert>
00016 
00017 class PixelClusterShapeSeedComparitor : public SeedComparitor {
00018     public:
00019         PixelClusterShapeSeedComparitor(const edm::ParameterSet &cfg) ;
00020         virtual ~PixelClusterShapeSeedComparitor() ; 
00021         virtual void init(const edm::EventSetup& es) ;
00022         virtual bool compatible(const SeedingHitSet  &hits, const TrackingRegion & region) const { return true; }
00023         virtual bool compatible(const TrajectorySeed &seed) const { return true; }
00024         virtual bool compatible(const TrajectoryStateOnSurface &,
00025                 const TransientTrackingRecHit::ConstRecHitPointer &hit) const ;
00026         virtual bool compatible(const SeedingHitSet  &hits, 
00027                 const GlobalTrajectoryParameters &helixStateAtVertex,
00028                 const FastHelix                  &helix,
00029                 const TrackingRegion & region) const ;
00030         virtual bool compatible(const SeedingHitSet  &hits, 
00031                 const GlobalTrajectoryParameters &straightLineStateAtVertex,
00032                 const TrackingRegion & region) const ;
00033 
00034     private:
00035         bool compatibleHit(const TransientTrackingRecHit &hit, const GlobalVector &direction) const ;
00036 
00037         std::string filterName_;
00038         mutable edm::ESHandle<ClusterShapeHitFilter> filterHandle_;
00039         bool filterAtHelixStage_;
00040         bool filterPixelHits_, filterStripHits_;
00041 };
00042 
00043 
00044 PixelClusterShapeSeedComparitor::PixelClusterShapeSeedComparitor(const edm::ParameterSet &cfg) :
00045     filterName_(cfg.getParameter<std::string>("ClusterShapeHitFilterName")),
00046     filterAtHelixStage_(cfg.getParameter<bool>("FilterAtHelixStage")),
00047     filterPixelHits_(cfg.getParameter<bool>("FilterPixelHits")),
00048     filterStripHits_(cfg.getParameter<bool>("FilterStripHits"))
00049 {
00050 }
00051 
00052 PixelClusterShapeSeedComparitor::~PixelClusterShapeSeedComparitor() 
00053 {
00054 }
00055 
00056 void
00057 PixelClusterShapeSeedComparitor::init(const edm::EventSetup& es) {
00058     es.get<CkfComponentsRecord>().get(filterName_, filterHandle_);
00059 }
00060 
00061 
00062 bool
00063 PixelClusterShapeSeedComparitor::compatible(const TrajectoryStateOnSurface &tsos,
00064                 const TransientTrackingRecHit::ConstRecHitPointer &hit) const
00065 {
00066     if (filterAtHelixStage_) return true;
00067     assert(hit->isValid() && tsos.isValid());
00068     return compatibleHit(*hit, tsos.globalDirection());
00069 }
00070 
00071 bool
00072 PixelClusterShapeSeedComparitor::compatible(const SeedingHitSet  &hits, 
00073         const GlobalTrajectoryParameters &straightLineStateAtVertex,
00074         const TrackingRegion & region) const 
00075 { 
00076     return true; 
00077 }
00078 
00079 bool
00080 PixelClusterShapeSeedComparitor::compatible(const SeedingHitSet  &hits, 
00081         const GlobalTrajectoryParameters &helixStateAtVertex,
00082         const FastHelix                  &helix,
00083         const TrackingRegion & region) const 
00084 { 
00085     if (!filterAtHelixStage_) return true;
00086     float xc = helix.circle().x0(), yc = helix.circle().y0();
00087 
00088     GlobalPoint  vertex = helixStateAtVertex.position();
00089     GlobalVector momvtx = helixStateAtVertex.momentum();
00090     float x0 = vertex.x(), y0 = vertex.y();
00091     for (unsigned int i = 0, n = hits.size(); i < n; ++i) {
00092         const TransientTrackingRecHit &hit = *hits[i];
00093         GlobalPoint pos = hit.globalPosition();
00094         float x1 = pos.x(), y1 = pos.y(), dx1 = x1 - xc, dy1 = y1 - yc;
00095 
00096         // now figure out the proper tangent vector
00097         float perpx = -dy1, perpy = dx1;
00098         if (perpx * (x1-x0) + perpy * (y1 - y0) < 0) {
00099             perpy = -perpy; perpx = -perpx;
00100         }
00101        
00102         // now normalize (perpx, perpy, 1.0) to momentum (px, py, pz)
00103         float perp2 = perpx*perpx + perpy*perpy; 
00104         float pmom2 = momvtx.x()*momvtx.x() + momvtx.y()*momvtx.y(), momz2 = momvtx.z()*momvtx.z(), mom2 = pmom2 + momz2;
00105         float perpscale = sqrt(pmom2/mom2 / perp2), zscale = sqrt((1-pmom2/mom2));
00106         GlobalVector gdir(perpx*perpscale, perpy*perpscale, (momvtx.z() > 0 ? zscale : -zscale));
00107 
00108         if (!compatibleHit(hit, gdir)) {
00109             return false; // not yet
00110         }
00111     }
00112     return true; 
00113 }
00114 
00115 bool 
00116 PixelClusterShapeSeedComparitor::compatibleHit(const TransientTrackingRecHit &hit, const GlobalVector &direction) const 
00117 {
00118     if (hit.geographicalId().subdetId() <= 2) {
00119         if (!filterPixelHits_) return true;    
00120         const SiPixelRecHit *pixhit = dynamic_cast<const SiPixelRecHit *>(hit.hit());
00121         if (pixhit == 0) throw cms::Exception("LogicError", "Found a valid hit on the pixel detector which is not a SiPixelRecHit\n");
00122         //printf("Cheching hi hit on detid %10d, local direction is x = %9.6f, y = %9.6f, z = %9.6f\n", hit.geographicalId().rawId(), direction.x(), direction.y(), direction.z());
00123         return filterHandle_->isCompatible(*pixhit, direction);
00124     } else {
00125         if (!filterStripHits_) return true;
00126         const std::type_info &tid = typeid(*hit.hit());
00127         if (tid == typeid(SiStripMatchedRecHit2D)) {
00128             const SiStripMatchedRecHit2D* matchedHit = dynamic_cast<const SiStripMatchedRecHit2D *>(hit.hit());
00129             assert(matchedHit != 0);
00130             return (filterHandle_->isCompatible(DetId(matchedHit->monoId()), matchedHit->monoCluster(), direction) &&
00131                     filterHandle_->isCompatible(DetId(matchedHit->stereoId()), matchedHit->stereoCluster(), direction));
00132         } else if (tid == typeid(SiStripRecHit2D)) {
00133             const SiStripRecHit2D* recHit = dynamic_cast<const SiStripRecHit2D *>(hit.hit());
00134             assert(recHit != 0);
00135             return filterHandle_->isCompatible(*recHit, direction);
00136         } else if (tid == typeid(ProjectedSiStripRecHit2D)) {
00137             const ProjectedSiStripRecHit2D* precHit = dynamic_cast<const ProjectedSiStripRecHit2D *>(hit.hit());
00138             assert(precHit != 0);
00139             return filterHandle_->isCompatible(precHit->originalHit(), direction);
00140         } else {
00141             //printf("Questo e' un %s, che ci fo?\n", tid.name());
00142             return true;
00143         }
00144     }
00145 }
00146 
00147 DEFINE_EDM_PLUGIN(SeedComparitorFactory, PixelClusterShapeSeedComparitor, "PixelClusterShapeSeedComparitor");