CMS 3D CMS Logo

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

Go to the documentation of this file.
00001 
00009 #include "HLTrigger/Egamma/interface/HLTElectronOneOEMinusOneOPFilterRegional.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 HLTElectronOneOEMinusOneOPFilterRegional::HLTElectronOneOEMinusOneOPFilterRegional(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    barrelcut_  = iConfig.getParameter<double> ("barrelcut");
00034    endcapcut_  = iConfig.getParameter<double> ("endcapcut");
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 
00047 HLTElectronOneOEMinusOneOPFilterRegional::~HLTElectronOneOEMinusOneOPFilterRegional(){}
00048 
00049 
00050 // ------------ method called to produce the data  ------------
00051 bool
00052 HLTElectronOneOEMinusOneOPFilterRegional::filter(edm::Event& iEvent, const edm::EventSetup& iSetup)
00053 {
00054 
00055   // The filter object
00056   using namespace trigger;
00057     std::auto_ptr<trigger::TriggerFilterObjectWithRefs> filterproduct (new trigger::TriggerFilterObjectWithRefs(path(),module()));
00058     if( store_ ){filterproduct->addCollectionTag(electronIsolatedProducer_);}
00059     if( store_ && !doIsolated_){filterproduct->addCollectionTag(electronNonIsolatedProducer_);}  
00060     //will be a collection of Ref<reco::ElectronCollection> ref;
00061     
00062   edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
00063   iEvent.getByLabel (candTag_,PrevFilterOutput);
00064 
00065   std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > recoecalcands;
00066   PrevFilterOutput->getObjects(TriggerCluster, recoecalcands);
00067   if(recoecalcands.empty()) PrevFilterOutput->getObjects(TriggerPhoton,recoecalcands);
00068 
00069   edm::Handle<reco::ElectronCollection> electronIsolatedHandle;
00070   iEvent.getByLabel(electronIsolatedProducer_,electronIsolatedHandle);
00071 
00072   edm::Handle<reco::ElectronCollection> electronNonIsolatedHandle;
00073   if(!doIsolated_) {
00074     iEvent.getByLabel(electronNonIsolatedProducer_,electronNonIsolatedHandle);
00075   }
00076 
00077  // look at all candidates,  check cuts and add to filter object
00078   int n(0);
00079   
00080     //loop over all the RecoCandidates from the previous filter, 
00081     // associate them with the corresponding Electron object 
00082     //(the matching is done checking the super clusters)
00083     // and put into the event a Ref to the Electron objects that passes the
00084     // selections  
00085  
00086   edm::RefToBase<reco::Candidate> candref;   
00087   for (unsigned int i=0; i<recoecalcands.size(); i++) {
00088     
00089     reco::SuperClusterRef recr2 = recoecalcands[i]->superCluster();
00090     //loop over the electrons to find the matching one
00091     for(reco::ElectronCollection::const_iterator iElectron = electronIsolatedHandle->begin(); iElectron != electronIsolatedHandle->end(); iElectron++){
00092       // ElectronRef is a Ref<reco::RecoEcalCandidateCollection>   
00093       reco::ElectronRef electronref(reco::ElectronRef(electronIsolatedHandle,iElectron - electronIsolatedHandle->begin()));
00094       const reco::SuperClusterRef theClus = electronref->superCluster();
00095       if(&(*recr2) ==  &(*theClus)) {
00096         
00097         float elecEoverp = 0;
00098         const math::XYZVector trackMom =  electronref->track()->momentum();
00099         if( trackMom.R() != 0) elecEoverp = 
00100                                  fabs((1/electronref->superCluster()->energy()) - (1/trackMom.R()));
00101 
00102         if( fabs(electronref->eta()) < 1.5 ){
00103           if ( elecEoverp < barrelcut_) {
00104             n++;
00105             filterproduct->addObject(TriggerElectron, electronref);
00106           }
00107         }
00108         if( fabs(electronref->eta()) > 1.5 ) {
00109           if ( elecEoverp < endcapcut_) {
00110             n++;
00111             filterproduct->addObject(TriggerElectron, electronref);
00112           }
00113         }
00114       }//end of the if checking the matching of the SC from RecoCandidate and the one from Electrons
00115     }//end of loop over electrons
00116 
00117     if(!doIsolated_) {
00118     //loop over the electrons to find the matching one
00119     for(reco::ElectronCollection::const_iterator iElectron = electronNonIsolatedHandle->begin(); iElectron != electronNonIsolatedHandle->end(); iElectron++){
00120       
00121       reco::ElectronRef electronref(reco::ElectronRef(electronNonIsolatedHandle,iElectron - electronNonIsolatedHandle->begin()));
00122       const reco::SuperClusterRef theClus = electronref->superCluster();
00123       
00124       if(&(*recr2) ==  &(*theClus)) {
00125 
00126         float elecEoverp = 0;
00127         const math::XYZVector trackMom =  electronref->track()->momentum();
00128         if( trackMom.R() != 0) elecEoverp = 
00129                                 fabs((1/electronref->superCluster()->energy()) - (1/trackMom.R())); 
00130 
00131         if( fabs(electronref->eta()) < 1.5 ){
00132           if ( elecEoverp < barrelcut_) {
00133             n++;
00134             filterproduct->addObject(TriggerElectron, electronref);
00135           }
00136         }
00137         if( fabs(electronref->eta()) > 1.5 ){
00138           if ( elecEoverp < endcapcut_) {
00139             n++;
00140             filterproduct->addObject(TriggerElectron, electronref);
00141           }
00142         }
00143       }//end of the if checking the matching of the SC from RecoCandidate and the one from Electrons
00144     }//end of loop over electrons
00145     }
00146   }//end of loop ober candidates
00147 
00148   // filter decision
00149   bool accept(n>=ncandcut_);
00150   
00151   // put filter object into the Event
00152   iEvent.put(filterproduct);
00153   
00154   return accept;
00155 }