CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/DPGAnalysis/Skims/src/FilterScrapingPixelProbability.cc

Go to the documentation of this file.
00001 
00002 // Original Author: Gavril Giurgiu (JHU) 
00003 
00004 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
00005 #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
00006 #include "DataFormats/SiPixelDetId/interface/PXFDetId.h"
00007 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00008 
00009 #include "DPGAnalysis/Skims/interface/FilterScrapingPixelProbability.h"
00010 
00011 using namespace edm;
00012 using namespace std;
00013 
00014 FilterScrapingPixelProbability::FilterScrapingPixelProbability(const edm::ParameterSet& iConfig)
00015 {
00016   //std::cout << "FilterScrapingPixelProbability::FilterScrapingPixelProbability " << std::endl;
00017 
00018   apply_filter                 = iConfig.getUntrackedParameter<bool>  ( "apply_filter"                , true  );
00019   select_collision             = iConfig.getUntrackedParameter<bool>  ( "select_collision"            , true  );
00020   select_pkam                  = iConfig.getUntrackedParameter<bool>  ( "select_pkam"                 , false );
00021   select_other                 = iConfig.getUntrackedParameter<bool>  ( "select_other"                , false );
00022   low_probability              = iConfig.getUntrackedParameter<double>( "low_probability"             , 0.0   );
00023   low_probability_fraction_cut = iConfig.getUntrackedParameter<double>( "low_probability_fraction_cut", 0.4   ); 
00024   tracks_ = iConfig.getUntrackedParameter<edm::InputTag>("src",edm::InputTag("generalTracks"));
00025 }
00026 
00027 FilterScrapingPixelProbability::~FilterScrapingPixelProbability()
00028 {
00029 }
00030 
00031 bool FilterScrapingPixelProbability::filter( edm::Event& iEvent, const edm::EventSetup& iSetup)
00032 {
00033   bool accepted = false;
00034  
00035   double low_probability_fraction = -999999.9;  
00036  
00037   float n_hits_low_prob = 0.0;
00038   float n_hits_barrel   = 0.0;
00039 
00040   // Loop over generalTracks
00041   edm::Handle<reco::TrackCollection> trackCollection;
00042   iEvent.getByLabel(tracks_, trackCollection);
00043   const reco::TrackCollection *tracks = trackCollection.product();
00044   reco::TrackCollection::const_iterator tciter;
00045 
00046   //std::cout << "(int)tracks->size() = " << (int)tracks->size() << std::endl;
00047   
00048   if ( (int)tracks->size() > 0 )
00049     {
00050       // Loop on tracks
00051       for ( tciter=tracks->begin(); tciter!=tracks->end(); ++tciter)
00052         {
00053           // First loop on hits: find pixel hits
00054           for ( trackingRecHit_iterator it = tciter->recHitsBegin(); it != tciter->recHitsEnd(); ++it) 
00055             {
00056               const TrackingRecHit &thit = **it;
00057               
00058               const SiPixelRecHit* matchedhit = dynamic_cast<const SiPixelRecHit*>(&thit);
00059               
00060               // Check if the RecHit is a PixelRecHit
00061               if ( matchedhit ) 
00062                 {
00063                   DetId detId = (*it)->geographicalId();
00064                   
00065                   // Only consider barrel pixel hits
00066                   if ( (int)detId.subdetId() == (int)PixelSubdetector::PixelBarrel ) 
00067                     {
00068                       n_hits_barrel = n_hits_barrel + 1.0;
00069 
00070                       double pixel_hit_probability = matchedhit->clusterProbability(0);
00071                        
00072                       if ( pixel_hit_probability <= 0.0 )
00073                         n_hits_low_prob = n_hits_low_prob + 1.0; 
00074                     
00075                     } // if ( (int)detId.subdetId() == (int)PixelSubdetector::PixelBarrel ) 
00076 
00077                 } // if ( matchedhit )
00078 
00079             } // for ( trackingRecHit_iterator it = tciter->recHitsBegin(); it != tciter->recHitsEnd(); ++it) 
00080           
00081         } //  for ( tciter=tracks->begin(); tciter!=tracks->end(); ++tciter) 
00082       
00083     } // if ( (int)tracks->size() > 0 )
00084 
00085   bool is_collision = false;
00086   bool is_pkam      = false;
00087   bool is_other     = false;
00088 
00089   if ( n_hits_barrel > 0.0 )
00090     {
00091       low_probability_fraction = n_hits_low_prob / n_hits_barrel;
00092 
00093       if ( low_probability_fraction < 0.4 )
00094         is_collision = true;
00095       else 
00096         is_pkam = true;
00097     }
00098   else
00099     is_other = true;
00100   
00101   if ( ( select_collision && is_collision ) || 
00102        ( select_pkam      && is_pkam      ) ||
00103        ( select_other     && is_other     ) )
00104     accepted = true;
00105  
00106   if ( apply_filter )
00107     return accepted;
00108   else
00109     return true;
00110 
00111 }
00112 
00113 //define this as a plug-in
00114 DEFINE_FWK_MODULE(FilterScrapingPixelProbability);