CMS 3D CMS Logo

FilterOutScraping.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: BeamSplash
4 // Class: BeamSPlash
5 //
6 //
7 // Original Author: Luca Malgeri
8 
9 #include <memory>
10 #include <vector>
11 #include <map>
12 #include <set>
13 
14 // user include files
15 
24 
26 
27 
28 using namespace edm;
29 using namespace std;
30 
32 {
33 
34  applyfilter = iConfig.getUntrackedParameter<bool>("applyfilter",true);
35  debugOn = iConfig.getUntrackedParameter<bool>("debugOn",false);
36  thresh = iConfig.getUntrackedParameter<double>("thresh",0.2);
37  numtrack = iConfig.getUntrackedParameter<unsigned int>("numtrack",10);
38  tracks_ = consumes<reco::TrackCollection>(iConfig.getUntrackedParameter<edm::InputTag>("src",edm::InputTag("generalTracks")));
39 }
40 
42 {
43 }
44 
46 {
47 
48  bool accepted = false;
49  float fraction = 0;
50  // get GeneralTracks collection
51 
53  iEvent.getByToken(tracks_,tkRef);
54  const reco::TrackCollection* tkColl = tkRef.product();
55 
56  //std::cout << "Total Number of Tracks " << tkColl->size() << std::endl;
57 
58  int numhighpurity=0;
59  _trackQuality = reco::TrackBase::qualityByName("highPurity");
60 
61  if(tkColl->size()>numtrack){
62  reco::TrackCollection::const_iterator itk = tkColl->begin();
63  reco::TrackCollection::const_iterator itk_e = tkColl->end();
64  for(;itk!=itk_e;++itk){
65  // std::cout << "HighPurity? " << itk->quality(_trackQuality) << std::endl;
66  if(itk->quality(_trackQuality)) numhighpurity++;
67  }
68  fraction = (float)numhighpurity/(float)tkColl->size();
69  if(fraction>thresh) accepted=true;
70  }else{
71  //if less than 10 Tracks accept the event anyway
72  accepted= true;
73  }
74 
75 
76  if (debugOn) {
77  int ievt = iEvent.id().event();
78  int irun = iEvent.id().run();
79  int ils = iEvent.luminosityBlock();
80  int bx = iEvent.bunchCrossing();
81 
82  std::cout << "FilterOutScraping_debug: Run " << irun << " Event " << ievt << " Lumi Block " << ils << " Bunch Crossing " << bx << " Fraction " << fraction << " NTracks " << tkColl->size() << " Accepted " << accepted << std::endl;
83  }
84 
85  if (applyfilter)
86  return accepted;
87  else
88  return true;
89 
90 }
91 
92 //define this as a plug-in
RunNumber_t run() const
Definition: EventID.h:39
EventNumber_t event() const
Definition: EventID.h:41
T getUntrackedParameter(std::string const &, T const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:15
int bunchCrossing() const
Definition: EventBase.h:64
edm::LuminosityBlockNumber_t luminosityBlock() const
Definition: EventBase.h:61
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:134
~FilterOutScraping() override
T const * product() const
Definition: Handle.h:74
bool accepted(std::vector< std::string_view > const &, std::string_view)
edm::EventID id() const
Definition: EventBase.h:59
HLT enums.
const std::vector< reco::CandidatePtr > & tracks_
bool filter(edm::Event &, const edm::EventSetup &) override
FilterOutScraping(const edm::ParameterSet &)