CMS 3D CMS Logo

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