CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
FilterScrapingPixelProbability Class Reference

#include <FilterScrapingPixelProbability.h>

Inheritance diagram for FilterScrapingPixelProbability:
edm::EDFilter edm::ProducerBase edm::ProductRegistryHelper

Public Member Functions

 FilterScrapingPixelProbability (const edm::ParameterSet &)
 
 ~FilterScrapingPixelProbability ()
 
- Public Member Functions inherited from edm::EDFilter
 EDFilter ()
 
virtual ~EDFilter ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
boost::function< void(const
BranchDescription &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 

Private Member Functions

virtual bool filter (edm::Event &, const edm::EventSetup &)
 

Private Attributes

bool apply_filter
 
double low_probability
 
double low_probability_fraction_cut
 
bool select_collision
 
bool select_other
 
bool select_pkam
 
edm::InputTag tracks_
 

Additional Inherited Members

- Public Types inherited from edm::EDFilter
typedef EDFilter ModuleType
 
typedef WorkerT< EDFilterWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDFilter
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDFilter
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
template<class TProducer , class TMethod >
void callWhenNewProductsRegistered (TProducer *iProd, TMethod iMethod)
 

Detailed Description

Definition at line 45 of file FilterScrapingPixelProbability.h.

Constructor & Destructor Documentation

FilterScrapingPixelProbability::FilterScrapingPixelProbability ( const edm::ParameterSet iConfig)
explicit

Definition at line 14 of file FilterScrapingPixelProbability.cc.

References skim_scrape_pixel_prob_cfg::apply_filter, funct::false, edm::ParameterSet::getUntrackedParameter(), skim_scrape_pixel_prob_cfg::low_probability, skim_scrape_pixel_prob_cfg::low_probability_fraction_cut, skim_scrape_pixel_prob_cfg::select_collision, skim_scrape_pixel_prob_cfg::select_other, skim_scrape_pixel_prob_cfg::select_pkam, tracks_, and funct::true.

15 {
16  //std::cout << "FilterScrapingPixelProbability::FilterScrapingPixelProbability " << std::endl;
17 
18  apply_filter = iConfig.getUntrackedParameter<bool> ( "apply_filter" , true );
19  select_collision = iConfig.getUntrackedParameter<bool> ( "select_collision" , true );
20  select_pkam = iConfig.getUntrackedParameter<bool> ( "select_pkam" , false );
21  select_other = iConfig.getUntrackedParameter<bool> ( "select_other" , false );
22  low_probability = iConfig.getUntrackedParameter<double>( "low_probability" , 0.0 );
23  low_probability_fraction_cut = iConfig.getUntrackedParameter<double>( "low_probability_fraction_cut", 0.4 );
24  tracks_ = iConfig.getUntrackedParameter<edm::InputTag>("src",edm::InputTag("generalTracks"));
25 }
T getUntrackedParameter(std::string const &, T const &) const
FilterScrapingPixelProbability::~FilterScrapingPixelProbability ( )

Definition at line 27 of file FilterScrapingPixelProbability.cc.

28 {
29 }

Member Function Documentation

bool FilterScrapingPixelProbability::filter ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
privatevirtual

Implements edm::EDFilter.

Definition at line 31 of file FilterScrapingPixelProbability.cc.

References skim_scrape_pixel_prob_cfg::apply_filter, SiPixelRecHit::clusterProbability(), edm::Event::getByLabel(), PixelSubdetector::PixelBarrel, edm::Handle< T >::product(), skim_scrape_pixel_prob_cfg::select_collision, skim_scrape_pixel_prob_cfg::select_other, skim_scrape_pixel_prob_cfg::select_pkam, DetId::subdetId(), testEve_cfg::tracks, and tracks_.

Referenced by Vispa.Plugins.Browser.BrowserTabController.BrowserTabController::filter(), Vispa.Plugins.Browser.BrowserTabController.BrowserTabController::find(), Vispa.Plugins.Browser.BrowserTabController.BrowserTabController::setDataAccessor(), and Vispa.Plugins.Browser.BrowserTabController.BrowserTabController::switchCenterView().

32 {
33  bool accepted = false;
34 
35  double low_probability_fraction = -999999.9;
36 
37  float n_hits_low_prob = 0.0;
38  float n_hits_barrel = 0.0;
39 
40  // Loop over generalTracks
41  edm::Handle<reco::TrackCollection> trackCollection;
42  iEvent.getByLabel(tracks_, trackCollection);
43  const reco::TrackCollection *tracks = trackCollection.product();
44  reco::TrackCollection::const_iterator tciter;
45 
46  //std::cout << "(int)tracks->size() = " << (int)tracks->size() << std::endl;
47 
48  if ( (int)tracks->size() > 0 )
49  {
50  // Loop on tracks
51  for ( tciter=tracks->begin(); tciter!=tracks->end(); ++tciter)
52  {
53  // First loop on hits: find pixel hits
54  for ( trackingRecHit_iterator it = tciter->recHitsBegin(); it != tciter->recHitsEnd(); ++it)
55  {
56  const TrackingRecHit &thit = **it;
57 
58  const SiPixelRecHit* matchedhit = dynamic_cast<const SiPixelRecHit*>(&thit);
59 
60  // Check if the RecHit is a PixelRecHit
61  if ( matchedhit )
62  {
63  DetId detId = (*it)->geographicalId();
64 
65  // Only consider barrel pixel hits
66  if ( (int)detId.subdetId() == (int)PixelSubdetector::PixelBarrel )
67  {
68  n_hits_barrel = n_hits_barrel + 1.0;
69 
70  double pixel_hit_probability = matchedhit->clusterProbability(0);
71 
72  if ( pixel_hit_probability <= 0.0 )
73  n_hits_low_prob = n_hits_low_prob + 1.0;
74 
75  } // if ( (int)detId.subdetId() == (int)PixelSubdetector::PixelBarrel )
76 
77  } // if ( matchedhit )
78 
79  } // for ( trackingRecHit_iterator it = tciter->recHitsBegin(); it != tciter->recHitsEnd(); ++it)
80 
81  } // for ( tciter=tracks->begin(); tciter!=tracks->end(); ++tciter)
82 
83  } // if ( (int)tracks->size() > 0 )
84 
85  bool is_collision = false;
86  bool is_pkam = false;
87  bool is_other = false;
88 
89  if ( n_hits_barrel > 0.0 )
90  {
91  low_probability_fraction = n_hits_low_prob / n_hits_barrel;
92 
93  if ( low_probability_fraction < 0.4 )
94  is_collision = true;
95  else
96  is_pkam = true;
97  }
98  else
99  is_other = true;
100 
101  if ( ( select_collision && is_collision ) ||
102  ( select_pkam && is_pkam ) ||
103  ( select_other && is_other ) )
104  accepted = true;
105 
106  if ( apply_filter )
107  return accepted;
108  else
109  return true;
110 
111 }
float clusterProbability(unsigned int flags=0) const
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:10
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:39
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
Definition: DetId.h:20
tuple tracks
Definition: testEve_cfg.py:39
T const * product() const
Definition: Handle.h:74
Our base class.
Definition: SiPixelRecHit.h:22

Member Data Documentation

bool FilterScrapingPixelProbability::apply_filter
private

Definition at line 54 of file FilterScrapingPixelProbability.h.

double FilterScrapingPixelProbability::low_probability
private

Definition at line 58 of file FilterScrapingPixelProbability.h.

double FilterScrapingPixelProbability::low_probability_fraction_cut
private

Definition at line 59 of file FilterScrapingPixelProbability.h.

bool FilterScrapingPixelProbability::select_collision
private

Definition at line 55 of file FilterScrapingPixelProbability.h.

bool FilterScrapingPixelProbability::select_other
private

Definition at line 57 of file FilterScrapingPixelProbability.h.

bool FilterScrapingPixelProbability::select_pkam
private

Definition at line 56 of file FilterScrapingPixelProbability.h.

edm::InputTag FilterScrapingPixelProbability::tracks_
private

Definition at line 61 of file FilterScrapingPixelProbability.h.