CMS 3D CMS Logo

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

Go to the documentation of this file.
00001 
00009 #include "HLTrigger/Egamma/interface/HLTElectronEoverpFilterRegional.h"
00010 
00011 #include "DataFormats/Common/interface/Handle.h"
00012 
00013 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
00014 
00015 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
00016 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidateFwd.h"
00017 
00018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00019 
00020 #include "DataFormats/EgammaCandidates/interface/Electron.h"
00021 #include "DataFormats/EgammaCandidates/interface/ElectronFwd.h"
00022 #include "DataFormats/TrackReco/interface/Track.h"
00023 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00024 
00025 //
00026 // constructors and destructor
00027 //
00028 HLTElectronEoverpFilterRegional::HLTElectronEoverpFilterRegional(const edm::ParameterSet& iConfig)
00029 {
00030    candTag_ = iConfig.getParameter< edm::InputTag > ("candTag");
00031    electronIsolatedProducer_ = iConfig.getParameter< edm::InputTag > ("electronIsolatedProducer");
00032    electronNonIsolatedProducer_ = iConfig.getParameter< edm::InputTag > ("electronNonIsolatedProducer");
00033    eoverpbarrelcut_  = iConfig.getParameter<double> ("eoverpbarrelcut");
00034    eoverpendcapcut_  = iConfig.getParameter<double> ("eoverpendcapcut");
00035    ncandcut_  = iConfig.getParameter<int> ("ncandcut");
00036    doIsolated_  = iConfig.getParameter<bool> ("doIsolated");
00037    
00038    store_ = iConfig.getUntrackedParameter<bool> ("SaveTag",false) ;
00039    // L1IsoCollTag_= iConfig.getParameter< edm::InputTag > ("L1IsoCand"); 
00040    // L1NonIsoCollTag_= iConfig.getParameter< edm::InputTag > ("L1NonIsoCand"); 
00041 
00042    //register your products
00043    produces<trigger::TriggerFilterObjectWithRefs>();
00044 }
00045 
00046 HLTElectronEoverpFilterRegional::~HLTElectronEoverpFilterRegional(){}
00047 
00048 
00049 // ------------ method called to produce the data  ------------
00050 bool
00051 HLTElectronEoverpFilterRegional::filter(edm::Event& iEvent, const edm::EventSetup& iSetup)
00052 {
00053 
00054   // The filter object
00055   using namespace trigger;
00056     std::auto_ptr<trigger::TriggerFilterObjectWithRefs> filterproduct (new trigger::TriggerFilterObjectWithRefs(path(),module()));
00057     if( store_ ){filterproduct->addCollectionTag(electronIsolatedProducer_);}
00058     if( store_ && !doIsolated_){filterproduct->addCollectionTag(electronNonIsolatedProducer_);}   
00059     //will be a collection of Ref<reco::ElectronCollection> ref;
00060 
00061   edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
00062   iEvent.getByLabel (candTag_,PrevFilterOutput);
00063 
00064   std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > recoecalcands;
00065   PrevFilterOutput->getObjects(TriggerCluster, recoecalcands);
00066 
00067    // Get the HLT electrons from EgammaHLTPixelMatchElectronProducers
00068   edm::Handle<reco::ElectronCollection> electronIsolatedHandle;
00069   iEvent.getByLabel(electronIsolatedProducer_,electronIsolatedHandle);
00070 
00071   edm::Handle<reco::ElectronCollection> electronNonIsolatedHandle;
00072   if(!doIsolated_) {
00073     iEvent.getByLabel(electronNonIsolatedProducer_,electronNonIsolatedHandle);
00074   }
00075 
00076  // look at all candidates,  check cuts and add to filter object
00077   int n(0);
00078   
00079     //loop over all the RecoCandidates from the previous filter, 
00080     // associate them with the corresponding Electron object 
00081     //(the matching is done checking the super clusters)
00082     // and put into the event a Ref to the Electron objects that passes the
00083     // selections  
00084   for (unsigned int i=0; i<recoecalcands.size(); i++) {
00085     reco::SuperClusterRef recr2 = recoecalcands[i]->superCluster();
00086 
00087     //loop over the electrons to find the matching one
00088     for(reco::ElectronCollection::const_iterator iElectron = electronIsolatedHandle->begin(); iElectron != electronIsolatedHandle->end(); iElectron++){
00089       
00090       reco::ElectronRef electronref(reco::ElectronRef(electronIsolatedHandle,iElectron - electronIsolatedHandle->begin()));
00091       const reco::SuperClusterRef theClus = electronref->superCluster();
00092       
00093       if(&(*recr2) ==  &(*theClus)) {
00094 
00095         float elecEoverp = 0;
00096         const math::XYZVector trackMom =  electronref->track()->momentum();
00097         if( trackMom.R() != 0) elecEoverp = 
00098                                  electronref->superCluster()->energy()/ trackMom.R();
00099 
00100         if( fabs(electronref->eta()) < 1.5 ){
00101           if ( elecEoverp < eoverpbarrelcut_) {
00102             n++;
00103             filterproduct->addObject(TriggerElectron, electronref);
00104           }
00105         }
00106         if( fabs(electronref->eta()) > 1.5 ){
00107           if ( elecEoverp < eoverpendcapcut_) {
00108             n++;
00109             filterproduct->addObject(TriggerElectron, electronref);
00110           }
00111         }
00112       }//end of the if checking the matching of the SC from RecoCandidate and the one from Electrons
00113     }//end of loop over electrons
00114 
00115     if(!doIsolated_) {
00116     //loop over the electrons to find the matching one
00117     for(reco::ElectronCollection::const_iterator iElectron = electronNonIsolatedHandle->begin(); iElectron != electronNonIsolatedHandle->end(); iElectron++){
00118       
00119       reco::ElectronRef electronref(reco::ElectronRef(electronNonIsolatedHandle,iElectron - electronNonIsolatedHandle->begin()));
00120       const reco::SuperClusterRef theClus = electronref->superCluster();
00121       
00122       if(&(*recr2) ==  &(*theClus)) {
00123 
00124         float elecEoverp = 0;
00125         const math::XYZVector trackMom =  electronref->track()->momentum();
00126         if( trackMom.R() != 0) elecEoverp = 
00127                                  electronref->superCluster()->energy()/ trackMom.R();
00128 
00129         if( fabs(electronref->eta()) < 1.5 ){
00130           if ( elecEoverp < eoverpbarrelcut_) {
00131             n++;
00132             filterproduct->addObject(TriggerElectron, electronref);
00133           }
00134         }
00135         if( fabs(electronref->eta()) > 1.5 ){
00136           if ( elecEoverp < eoverpendcapcut_) {
00137             n++;
00138             filterproduct->addObject(TriggerElectron, electronref);
00139           }
00140         }
00141       }//end of the if checking the matching of the SC from RecoCandidate and the one from Electrons
00142     }//end of loop over electrons
00143     }
00144   }//end of loop ober candidates
00145 
00146   // filter decision
00147   bool accept(n>=ncandcut_);
00148   
00149   // put filter object into the Event
00150   iEvent.put(filterproduct);
00151   
00152   return accept;
00153 }