CMS 3D CMS Logo

Public Member Functions | Private Attributes

HLTElectronEoverpFilterRegional Class Reference

#include <HLTElectronEoverpFilterRegional.h>

Inheritance diagram for HLTElectronEoverpFilterRegional:
HLTFilter edm::EDFilter edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

virtual bool filter (edm::Event &, const edm::EventSetup &)
 HLTElectronEoverpFilterRegional (const edm::ParameterSet &)
 ~HLTElectronEoverpFilterRegional ()

Private Attributes

edm::InputTag candTag_
bool doIsolated_
edm::InputTag electronIsolatedProducer_
edm::InputTag electronNonIsolatedProducer_
double eoverpbarrelcut_
double eoverpendcapcut_
int ncandcut_
bool store_

Detailed Description

Author:
Monica Vazquez Acosta (CERN)
Id:
HLTElectronEoverpFilterRegional.h,v 1.3 2009/02/09 16:27:17 covarell Exp
Author:
Monica Vazquez Acosta (CERN)
Id:
HLTElectronEoverpFilterRegional.cc,v 1.9 2009/02/09 16:27:18 covarell Exp

Definition at line 16 of file HLTElectronEoverpFilterRegional.h.


Constructor & Destructor Documentation

HLTElectronEoverpFilterRegional::HLTElectronEoverpFilterRegional ( const edm::ParameterSet iConfig) [explicit]

Definition at line 28 of file HLTElectronEoverpFilterRegional.cc.

References candTag_, doIsolated_, electronIsolatedProducer_, electronNonIsolatedProducer_, eoverpbarrelcut_, eoverpendcapcut_, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), ncandcut_, and store_.

{
   candTag_ = iConfig.getParameter< edm::InputTag > ("candTag");
   electronIsolatedProducer_ = iConfig.getParameter< edm::InputTag > ("electronIsolatedProducer");
   electronNonIsolatedProducer_ = iConfig.getParameter< edm::InputTag > ("electronNonIsolatedProducer");
   eoverpbarrelcut_  = iConfig.getParameter<double> ("eoverpbarrelcut");
   eoverpendcapcut_  = iConfig.getParameter<double> ("eoverpendcapcut");
   ncandcut_  = iConfig.getParameter<int> ("ncandcut");
   doIsolated_  = iConfig.getParameter<bool> ("doIsolated");
   
   store_ = iConfig.getUntrackedParameter<bool> ("SaveTag",false) ;
   // L1IsoCollTag_= iConfig.getParameter< edm::InputTag > ("L1IsoCand"); 
   // L1NonIsoCollTag_= iConfig.getParameter< edm::InputTag > ("L1NonIsoCand"); 

   //register your products
   produces<trigger::TriggerFilterObjectWithRefs>();
}
HLTElectronEoverpFilterRegional::~HLTElectronEoverpFilterRegional ( )

Definition at line 46 of file HLTElectronEoverpFilterRegional.cc.

{}

Member Function Documentation

bool HLTElectronEoverpFilterRegional::filter ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [virtual]

Implements HLTFilter.

Definition at line 51 of file HLTElectronEoverpFilterRegional.cc.

References accept(), candTag_, doIsolated_, electronIsolatedProducer_, electronNonIsolatedProducer_, eoverpbarrelcut_, eoverpendcapcut_, edm::Event::getByLabel(), i, module(), n, ncandcut_, path(), edm::Event::put(), store_, trigger::TriggerCluster, and trigger::TriggerElectron.

{

  // The filter object
  using namespace trigger;
    std::auto_ptr<trigger::TriggerFilterObjectWithRefs> filterproduct (new trigger::TriggerFilterObjectWithRefs(path(),module()));
    if( store_ ){filterproduct->addCollectionTag(electronIsolatedProducer_);}
    if( store_ && !doIsolated_){filterproduct->addCollectionTag(electronNonIsolatedProducer_);}   
    //will be a collection of Ref<reco::ElectronCollection> ref;

  edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
  iEvent.getByLabel (candTag_,PrevFilterOutput);

  std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > recoecalcands;
  PrevFilterOutput->getObjects(TriggerCluster, recoecalcands);

   // Get the HLT electrons from EgammaHLTPixelMatchElectronProducers
  edm::Handle<reco::ElectronCollection> electronIsolatedHandle;
  iEvent.getByLabel(electronIsolatedProducer_,electronIsolatedHandle);

  edm::Handle<reco::ElectronCollection> electronNonIsolatedHandle;
  if(!doIsolated_) {
    iEvent.getByLabel(electronNonIsolatedProducer_,electronNonIsolatedHandle);
  }

 // look at all candidates,  check cuts and add to filter object
  int n(0);
  
    //loop over all the RecoCandidates from the previous filter, 
    // associate them with the corresponding Electron object 
    //(the matching is done checking the super clusters)
    // and put into the event a Ref to the Electron objects that passes the
    // selections  
  for (unsigned int i=0; i<recoecalcands.size(); i++) {
    reco::SuperClusterRef recr2 = recoecalcands[i]->superCluster();

    //loop over the electrons to find the matching one
    for(reco::ElectronCollection::const_iterator iElectron = electronIsolatedHandle->begin(); iElectron != electronIsolatedHandle->end(); iElectron++){
      
      reco::ElectronRef electronref(reco::ElectronRef(electronIsolatedHandle,iElectron - electronIsolatedHandle->begin()));
      const reco::SuperClusterRef theClus = electronref->superCluster();
      
      if(&(*recr2) ==  &(*theClus)) {

        float elecEoverp = 0;
        const math::XYZVector trackMom =  electronref->track()->momentum();
        if( trackMom.R() != 0) elecEoverp = 
                                 electronref->superCluster()->energy()/ trackMom.R();

        if( fabs(electronref->eta()) < 1.5 ){
          if ( elecEoverp < eoverpbarrelcut_) {
            n++;
            filterproduct->addObject(TriggerElectron, electronref);
          }
        }
        if( fabs(electronref->eta()) > 1.5 ){
          if ( elecEoverp < eoverpendcapcut_) {
            n++;
            filterproduct->addObject(TriggerElectron, electronref);
          }
        }
      }//end of the if checking the matching of the SC from RecoCandidate and the one from Electrons
    }//end of loop over electrons

    if(!doIsolated_) {
    //loop over the electrons to find the matching one
    for(reco::ElectronCollection::const_iterator iElectron = electronNonIsolatedHandle->begin(); iElectron != electronNonIsolatedHandle->end(); iElectron++){
      
      reco::ElectronRef electronref(reco::ElectronRef(electronNonIsolatedHandle,iElectron - electronNonIsolatedHandle->begin()));
      const reco::SuperClusterRef theClus = electronref->superCluster();
      
      if(&(*recr2) ==  &(*theClus)) {

        float elecEoverp = 0;
        const math::XYZVector trackMom =  electronref->track()->momentum();
        if( trackMom.R() != 0) elecEoverp = 
                                 electronref->superCluster()->energy()/ trackMom.R();

        if( fabs(electronref->eta()) < 1.5 ){
          if ( elecEoverp < eoverpbarrelcut_) {
            n++;
            filterproduct->addObject(TriggerElectron, electronref);
          }
        }
        if( fabs(electronref->eta()) > 1.5 ){
          if ( elecEoverp < eoverpendcapcut_) {
            n++;
            filterproduct->addObject(TriggerElectron, electronref);
          }
        }
      }//end of the if checking the matching of the SC from RecoCandidate and the one from Electrons
    }//end of loop over electrons
    }
  }//end of loop ober candidates

  // filter decision
  bool accept(n>=ncandcut_);
  
  // put filter object into the Event
  iEvent.put(filterproduct);
  
  return accept;
}

Member Data Documentation

Definition at line 24 of file HLTElectronEoverpFilterRegional.h.

Referenced by filter(), and HLTElectronEoverpFilterRegional().

Definition at line 27 of file HLTElectronEoverpFilterRegional.h.

Referenced by filter(), and HLTElectronEoverpFilterRegional().

Definition at line 25 of file HLTElectronEoverpFilterRegional.h.

Referenced by filter(), and HLTElectronEoverpFilterRegional().

Definition at line 26 of file HLTElectronEoverpFilterRegional.h.

Referenced by filter(), and HLTElectronEoverpFilterRegional().

Definition at line 28 of file HLTElectronEoverpFilterRegional.h.

Referenced by filter(), and HLTElectronEoverpFilterRegional().

Definition at line 29 of file HLTElectronEoverpFilterRegional.h.

Referenced by filter(), and HLTElectronEoverpFilterRegional().

Definition at line 30 of file HLTElectronEoverpFilterRegional.h.

Referenced by filter(), and HLTElectronEoverpFilterRegional().

Definition at line 32 of file HLTElectronEoverpFilterRegional.h.

Referenced by filter(), and HLTElectronEoverpFilterRegional().