CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/HLTrigger/Egamma/src/HLTElectronGenericFilter.cc

Go to the documentation of this file.
00001 
00009 #include "HLTrigger/Egamma/interface/HLTElectronGenericFilter.h"
00010 
00011 #include "DataFormats/Common/interface/Handle.h"
00012 
00013 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
00014 
00015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00016 
00017 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
00018 #include "DataFormats/EgammaCandidates/interface/Electron.h"
00019 #include "DataFormats/EgammaCandidates/interface/ElectronIsolationAssociation.h"
00020 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00021 #include "DataFormats/EgammaReco/interface/SuperClusterFwd.h"
00022 
00023 #include "DataFormats/Common/interface/AssociationMap.h"
00024 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
00025 
00026 //
00027 // constructors and destructor
00028 //
00029 HLTElectronGenericFilter::HLTElectronGenericFilter(const edm::ParameterSet& iConfig){
00030   candTag_ = iConfig.getParameter< edm::InputTag > ("candTag");
00031   isoTag_ = iConfig.getParameter< edm::InputTag > ("isoTag");
00032   nonIsoTag_ = iConfig.getParameter< edm::InputTag > ("nonIsoTag");
00033 
00034   lessThan_ = iConfig.getParameter<bool> ("lessThan");
00035   thrRegularEB_ = iConfig.getParameter<double> ("thrRegularEB");
00036   thrRegularEE_ = iConfig.getParameter<double> ("thrRegularEE");
00037   thrOverPtEB_ = iConfig.getParameter<double> ("thrOverPtEB");
00038   thrOverPtEE_ = iConfig.getParameter<double> ("thrOverPtEE");
00039   thrTimesPtEB_ = iConfig.getParameter<double> ("thrTimesPtEB");
00040   thrTimesPtEE_ = iConfig.getParameter<double> ("thrTimesPtEE");
00041   
00042   ncandcut_  = iConfig.getParameter<int> ("ncandcut");
00043   doIsolated_ = iConfig.getParameter<bool> ("doIsolated");
00044 
00045   store_ = iConfig.getUntrackedParameter<bool> ("SaveTag",false) ;
00046   L1IsoCollTag_= iConfig.getParameter< edm::InputTag > ("L1IsoCand"); 
00047   L1NonIsoCollTag_= iConfig.getParameter< edm::InputTag > ("L1NonIsoCand"); 
00048 
00049 //register your products
00050 produces<trigger::TriggerFilterObjectWithRefs>();
00051 }
00052 
00053 HLTElectronGenericFilter::~HLTElectronGenericFilter(){}
00054 
00055 
00056 // ------------ method called to produce the data  ------------
00057 bool
00058 HLTElectronGenericFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup)
00059 {
00060   using namespace trigger;
00061   std::auto_ptr<trigger::TriggerFilterObjectWithRefs> filterproduct (new trigger::TriggerFilterObjectWithRefs(path(),module()));
00062   if( store_ ){filterproduct->addCollectionTag(L1IsoCollTag_);}
00063   if( store_ && !doIsolated_){filterproduct->addCollectionTag(L1NonIsoCollTag_);}
00064 
00065   // Ref to Candidate object to be recorded in filter object
00066   reco::ElectronRef ref;
00067 
00068   edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
00069 
00070   iEvent.getByLabel (candTag_,PrevFilterOutput);
00071 
00072   std::vector<edm::Ref<reco::ElectronCollection> > elecands;
00073   PrevFilterOutput->getObjects(TriggerElectron, elecands);
00074 
00075   
00076   //get hold of isolated association map
00077   edm::Handle<reco::ElectronIsolationMap> depMap;
00078   iEvent.getByLabel (isoTag_,depMap);
00079   
00080   //get hold of non-isolated association map
00081   edm::Handle<reco::ElectronIsolationMap> depNonIsoMap;
00082   if(!doIsolated_) iEvent.getByLabel (nonIsoTag_,depNonIsoMap);
00083   
00084   // look at all photons, check cuts and add to filter object
00085   int n = 0;
00086   
00087   for (unsigned int i=0; i<elecands.size(); i++) {
00088     
00089     ref = elecands[i];
00090     reco::ElectronIsolationMap::const_iterator mapi = (*depMap).find( ref );    
00091     if (mapi==(*depMap).end() && !doIsolated_) mapi = (*depNonIsoMap).find( ref ); 
00092    
00093     float vali = mapi->val;
00094     float Pt = ref->pt();
00095     float Eta = fabs(ref->eta());
00096 
00097     if ( lessThan_ ) {
00098       if ( (Eta < 1.479 && vali <= thrRegularEB_) || (Eta >= 1.479 && vali <= thrRegularEE_) ) {
00099         n++;
00100         filterproduct->addObject(TriggerElectron, ref);
00101         continue;
00102       }
00103       if (Pt > 0. && (thrOverPtEB_ > 0. || thrOverPtEE_ > 0. || thrTimesPtEB_ > 0. || thrTimesPtEE_ > 0.) ) {
00104         if ((Eta < 1.479 && vali/Pt <= thrOverPtEB_) || (Eta >= 1.479 && vali/Pt <= thrOverPtEE_) ) {
00105           n++;
00106           filterproduct->addObject(TriggerElectron, ref);
00107           continue;
00108         }
00109         if ((Eta < 1.479 && vali*Pt <= thrTimesPtEB_) || (Eta >= 1.479 && vali*Pt <= thrTimesPtEE_) ) {
00110           n++;
00111           filterproduct->addObject(TriggerElectron, ref);
00112         }
00113       }
00114     } else {
00115       if ( (Eta < 1.479 && vali >= thrRegularEB_) || (Eta >= 1.479 && vali >= thrRegularEE_) ) {
00116         n++;
00117         filterproduct->addObject(TriggerElectron, ref);
00118         continue;
00119       }
00120       if (Pt > 0. && (thrOverPtEB_ > 0. || thrOverPtEE_ > 0. || thrTimesPtEB_ > 0. || thrTimesPtEE_ > 0.) ) {
00121         if ((Eta < 1.479 && vali/Pt >= thrOverPtEB_) || (Eta >= 1.479 && vali/Pt >= thrOverPtEE_) ) {
00122           n++;
00123           filterproduct->addObject(TriggerElectron, ref);
00124           continue;
00125         }
00126         if ((Eta < 1.479 && vali*Pt >= thrTimesPtEB_) || (Eta >= 1.479 && vali*Pt >= thrTimesPtEE_) ) {
00127           n++;
00128           filterproduct->addObject(TriggerElectron, ref);
00129         }
00130       }
00131     }
00132   }
00133   
00134   // filter decision
00135   bool accept(n>=ncandcut_);
00136   
00137   // put filter object into the Event
00138   iEvent.put(filterproduct);
00139 
00140   return accept;
00141 }
00142