File indexing completed on 2024-04-06 12:06:48
0001
0002
0003
0004 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
0005 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0006
0007 #include "DPGAnalysis/Skims/interface/FilterScrapingPixelProbability.h"
0008
0009 using namespace edm;
0010 using namespace std;
0011
0012 FilterScrapingPixelProbability::FilterScrapingPixelProbability(const edm::ParameterSet& iConfig) {
0013
0014
0015 apply_filter = iConfig.getUntrackedParameter<bool>("apply_filter", true);
0016 select_collision = iConfig.getUntrackedParameter<bool>("select_collision", true);
0017 select_pkam = iConfig.getUntrackedParameter<bool>("select_pkam", false);
0018 select_other = iConfig.getUntrackedParameter<bool>("select_other", false);
0019 low_probability = iConfig.getUntrackedParameter<double>("low_probability", 0.0);
0020 low_probability_fraction_cut = iConfig.getUntrackedParameter<double>("low_probability_fraction_cut", 0.4);
0021 tracks_ = iConfig.getUntrackedParameter<edm::InputTag>("src", edm::InputTag("generalTracks"));
0022 }
0023
0024 FilterScrapingPixelProbability::~FilterScrapingPixelProbability() {}
0025
0026 bool FilterScrapingPixelProbability::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0027 bool accepted = false;
0028
0029 double low_probability_fraction = -999999.9;
0030
0031 float n_hits_low_prob = 0.0;
0032 float n_hits_barrel = 0.0;
0033
0034
0035 edm::Handle<reco::TrackCollection> trackCollection;
0036 iEvent.getByLabel(tracks_, trackCollection);
0037 const reco::TrackCollection* tracks = trackCollection.product();
0038 reco::TrackCollection::const_iterator tciter;
0039
0040
0041
0042 if ((int)tracks->size() > 0) {
0043
0044 for (tciter = tracks->begin(); tciter != tracks->end(); ++tciter) {
0045
0046 for (trackingRecHit_iterator it = tciter->recHitsBegin(); it != tciter->recHitsEnd(); ++it) {
0047 const TrackingRecHit& thit = **it;
0048
0049 const SiPixelRecHit* matchedhit = dynamic_cast<const SiPixelRecHit*>(&thit);
0050
0051
0052 if (matchedhit) {
0053 DetId detId = (*it)->geographicalId();
0054
0055
0056 if ((int)detId.subdetId() == (int)PixelSubdetector::PixelBarrel) {
0057 n_hits_barrel = n_hits_barrel + 1.0;
0058
0059 double pixel_hit_probability = matchedhit->clusterProbability(0);
0060
0061 if (pixel_hit_probability <= 0.0)
0062 n_hits_low_prob = n_hits_low_prob + 1.0;
0063
0064 }
0065
0066 }
0067
0068 }
0069
0070 }
0071
0072 }
0073
0074 bool is_collision = false;
0075 bool is_pkam = false;
0076 bool is_other = false;
0077
0078 if (n_hits_barrel > 0.0) {
0079 low_probability_fraction = n_hits_low_prob / n_hits_barrel;
0080
0081 if (low_probability_fraction < 0.4)
0082 is_collision = true;
0083 else
0084 is_pkam = true;
0085 } else
0086 is_other = true;
0087
0088 if ((select_collision && is_collision) || (select_pkam && is_pkam) || (select_other && is_other))
0089 accepted = true;
0090
0091 if (apply_filter)
0092 return accepted;
0093 else
0094 return true;
0095 }
0096
0097
0098 DEFINE_FWK_MODULE(FilterScrapingPixelProbability);