CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FilterScrapingPixelProbability.cc
Go to the documentation of this file.
1 
2 // Original Author: Gavril Giurgiu (JHU)
3 
6 
8 
9 using namespace edm;
10 using namespace std;
11 
13 {
14  //std::cout << "FilterScrapingPixelProbability::FilterScrapingPixelProbability " << std::endl;
15 
16  apply_filter = iConfig.getUntrackedParameter<bool> ( "apply_filter" , true );
17  select_collision = iConfig.getUntrackedParameter<bool> ( "select_collision" , true );
18  select_pkam = iConfig.getUntrackedParameter<bool> ( "select_pkam" , false );
19  select_other = iConfig.getUntrackedParameter<bool> ( "select_other" , false );
20  low_probability = iConfig.getUntrackedParameter<double>( "low_probability" , 0.0 );
21  low_probability_fraction_cut = iConfig.getUntrackedParameter<double>( "low_probability_fraction_cut", 0.4 );
22  tracks_ = iConfig.getUntrackedParameter<edm::InputTag>("src",edm::InputTag("generalTracks"));
23 }
24 
26 {
27 }
28 
30 {
31  bool accepted = false;
32 
33  double low_probability_fraction = -999999.9;
34 
35  float n_hits_low_prob = 0.0;
36  float n_hits_barrel = 0.0;
37 
38  // Loop over generalTracks
40  iEvent.getByLabel(tracks_, trackCollection);
41  const reco::TrackCollection *tracks = trackCollection.product();
42  reco::TrackCollection::const_iterator tciter;
43 
44  //std::cout << "(int)tracks->size() = " << (int)tracks->size() << std::endl;
45 
46  if ( (int)tracks->size() > 0 )
47  {
48  // Loop on tracks
49  for ( tciter=tracks->begin(); tciter!=tracks->end(); ++tciter)
50  {
51  // First loop on hits: find pixel hits
52  for ( trackingRecHit_iterator it = tciter->recHitsBegin(); it != tciter->recHitsEnd(); ++it)
53  {
54  const TrackingRecHit &thit = **it;
55 
56  const SiPixelRecHit* matchedhit = dynamic_cast<const SiPixelRecHit*>(&thit);
57 
58  // Check if the RecHit is a PixelRecHit
59  if ( matchedhit )
60  {
61  DetId detId = (*it)->geographicalId();
62 
63  // Only consider barrel pixel hits
64  if ( (int)detId.subdetId() == (int)PixelSubdetector::PixelBarrel )
65  {
66  n_hits_barrel = n_hits_barrel + 1.0;
67 
68  double pixel_hit_probability = matchedhit->clusterProbability(0);
69 
70  if ( pixel_hit_probability <= 0.0 )
71  n_hits_low_prob = n_hits_low_prob + 1.0;
72 
73  } // if ( (int)detId.subdetId() == (int)PixelSubdetector::PixelBarrel )
74 
75  } // if ( matchedhit )
76 
77  } // for ( trackingRecHit_iterator it = tciter->recHitsBegin(); it != tciter->recHitsEnd(); ++it)
78 
79  } // for ( tciter=tracks->begin(); tciter!=tracks->end(); ++tciter)
80 
81  } // if ( (int)tracks->size() > 0 )
82 
83  bool is_collision = false;
84  bool is_pkam = false;
85  bool is_other = false;
86 
87  if ( n_hits_barrel > 0.0 )
88  {
89  low_probability_fraction = n_hits_low_prob / n_hits_barrel;
90 
91  if ( low_probability_fraction < 0.4 )
92  is_collision = true;
93  else
94  is_pkam = true;
95  }
96  else
97  is_other = true;
98 
99  if ( ( select_collision && is_collision ) ||
100  ( select_pkam && is_pkam ) ||
101  ( select_other && is_other ) )
102  accepted = true;
103 
104  if ( apply_filter )
105  return accepted;
106  else
107  return true;
108 
109 }
110 
111 //define this as a plug-in
T getUntrackedParameter(std::string const &, T const &) const
float clusterProbability(unsigned int flags=0) const
const std::vector< reco::PFCandidatePtr > & tracks_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
int iEvent
Definition: GenABIO.cc:230
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:420
Definition: DetId.h:18
T const * product() const
Definition: Handle.h:81
tuple tracks
Definition: testEve_cfg.py:39
volatile std::atomic< bool > shutdown_flag false
virtual bool filter(edm::Event &, const edm::EventSetup &) override
FilterScrapingPixelProbability(const edm::ParameterSet &)
TrackingRecHitCollection::base::const_iterator trackingRecHit_iterator
iterator over a vector of reference to TrackingRecHit in the same collection
Our base class.
Definition: SiPixelRecHit.h:23