CMS 3D CMS Logo

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

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:   BeamSplash
00004 // Class:     BeamSPlash
00005 //
00006 //
00007 // Original Author:  Luca Malgeri
00008 
00009 #include <memory>
00010 #include <vector>
00011 #include <map>
00012 #include <set>
00013 
00014 // user include files
00015 
00016 #include "FWCore/Utilities/interface/InputTag.h"
00017 #include "FWCore/Framework/interface/Frameworkfwd.h"
00018 #include "FWCore/Framework/interface/EDFilter.h"
00019 #include "FWCore/Framework/interface/Event.h"
00020 #include "FWCore/Framework/interface/MakerMacros.h"
00021 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00022 #include "FWCore/Framework/interface/ESHandle.h"
00023 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00024 
00025 #include "DPGAnalysis/Skims/interface/FilterOutScraping.h"
00026 
00027 
00028 using namespace edm;
00029 using namespace std;
00030 
00031 FilterOutScraping::FilterOutScraping(const edm::ParameterSet& iConfig)
00032 {
00033   
00034   applyfilter = iConfig.getUntrackedParameter<bool>("applyfilter",true);
00035   debugOn     = iConfig.getUntrackedParameter<bool>("debugOn",false);
00036   thresh =  iConfig.getUntrackedParameter<double>("thresh",0.2);
00037   numtrack = iConfig.getUntrackedParameter<unsigned int>("numtrack",10);
00038   tracks_ = iConfig.getUntrackedParameter<edm::InputTag>("src",edm::InputTag("generalTracks"));
00039 }
00040 
00041 FilterOutScraping::~FilterOutScraping()
00042 {
00043 }
00044 
00045 bool FilterOutScraping::filter( edm::Event& iEvent, const edm::EventSetup& iSetup)
00046 {
00047 
00048   bool accepted = false;
00049   float fraction = 0;  
00050   // get GeneralTracks collection
00051 
00052   edm::Handle<reco::TrackCollection> tkRef;
00053   iEvent.getByLabel(tracks_,tkRef);    
00054   const reco::TrackCollection* tkColl = tkRef.product();
00055 
00056   //std::cout << "Total Number of Tracks " << tkColl->size() << std::endl;
00057   
00058   int numhighpurity=0;
00059   _trackQuality = reco::TrackBase::qualityByName("highPurity");
00060 
00061   if(tkColl->size()>numtrack){ 
00062     reco::TrackCollection::const_iterator itk = tkColl->begin();
00063     reco::TrackCollection::const_iterator itk_e = tkColl->end();
00064     for(;itk!=itk_e;++itk){
00065       // std::cout << "HighPurity?  " << itk->quality(_trackQuality) << std::endl;
00066       if(itk->quality(_trackQuality)) numhighpurity++;
00067     }
00068     fraction = (float)numhighpurity/(float)tkColl->size();
00069     if(fraction>thresh) accepted=true;
00070   }else{
00071     //if less than 10 Tracks accept the event anyway    
00072     accepted= true;
00073   }
00074   
00075   
00076   if (debugOn) {
00077     int ievt = iEvent.id().event();
00078     int irun = iEvent.id().run();
00079     int ils = iEvent.luminosityBlock();
00080     int bx = iEvent.bunchCrossing();
00081     
00082     std::cout << "FilterOutScraping_debug: Run " << irun << " Event " << ievt << " Lumi Block " << ils << " Bunch Crossing " << bx << " Fraction " << fraction << " NTracks " << tkColl->size() << " Accepted " << accepted << std::endl;
00083   }
00084  
00085   if (applyfilter)
00086     return accepted;
00087   else
00088     return true;
00089 
00090 }
00091 
00092 //define this as a plug-in
00093 DEFINE_FWK_MODULE(FilterOutScraping);