CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

PFTauVertexSelector Class Reference

#include <PFTauVertexSelector.h>

Inheritance diagram for PFTauVertexSelector:
edm::EDFilter edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

 PFTauVertexSelector (const edm::ParameterSet &iConfig)
 ~PFTauVertexSelector ()

Private Member Functions

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

Private Attributes

edm::InputTag beamSpotSrc_
double dZ_
uint32_t filterOnNTaus_
std::vector< edm::InputTagrecoCandidateSrc_
edm::InputTag tauSrc_
std::vector< edm::InputTagtrackSrc_
edm::InputTag triggerFilterElectronsSrc_
edm::InputTag triggerFilterMuonsSrc_
bool useBeamSpot_
bool useLeadingRecoCandidate_
bool useLeadingTrack_
bool useTriggerFilterElectrons_
bool useTriggerFilterMuons_
bool useVertex_
edm::InputTag vertexSrc_

Detailed Description

Definition at line 23 of file PFTauVertexSelector.h.


Constructor & Destructor Documentation

PFTauVertexSelector::PFTauVertexSelector ( const edm::ParameterSet iConfig) [inline, explicit]

Definition at line 25 of file PFTauVertexSelector.h.

References beamSpotSrc_, dZ_, filterOnNTaus_, edm::ParameterSet::getParameter(), recoCandidateSrc_, tauSrc_, trackSrc_, triggerFilterElectronsSrc_, triggerFilterMuonsSrc_, useBeamSpot_, useLeadingRecoCandidate_, useLeadingTrack_, useTriggerFilterElectrons_, useTriggerFilterMuons_, useVertex_, and vertexSrc_.

                                                                  {   
         tauSrc_ = iConfig.getParameter<edm::InputTag>("tauSrc");
         useVertex_ = iConfig.getParameter<bool>("useVertex");
         vertexSrc_ = iConfig.getParameter<edm::InputTag>("vertexSrc");
         useBeamSpot_ = iConfig.getParameter<bool>("useBeamSpot");
         beamSpotSrc_ = iConfig.getParameter<edm::InputTag>("beamSpotSrc");
         useLeadingTrack_ = iConfig.getParameter<bool>("useLeadingTrack");
         trackSrc_ = iConfig.getParameter<std::vector<edm::InputTag> >("trackSrc");
         useLeadingRecoCandidate_ = iConfig.getParameter<bool>("useLeadingRecoCandidate");
         recoCandidateSrc_ = iConfig.getParameter<std::vector<edm::InputTag> >("recoCandidateSrc");
         useTriggerFilterElectrons_ = iConfig.getParameter<bool>("useTriggerFilterElectrons");
         triggerFilterElectronsSrc_ = iConfig.getParameter<edm::InputTag>("triggerFilterElectronsSrc");
         useTriggerFilterMuons_ = iConfig.getParameter<bool>("useTriggerFilterMuons");
         triggerFilterMuonsSrc_ = iConfig.getParameter<edm::InputTag>("triggerFilterMuonsSrc");
         dZ_ = iConfig.getParameter<double>("dZ");
         filterOnNTaus_ = iConfig.getParameter<uint32_t>("filterOnNTaus");
         produces<reco::PFTauCollection>();
      }
PFTauVertexSelector::~PFTauVertexSelector ( ) [inline]

Definition at line 43 of file PFTauVertexSelector.h.

{} 

Member Function Documentation

bool PFTauVertexSelector::filter ( edm::Event event,
const edm::EventSetup eventSetup 
) [private, virtual]

Implements edm::EDFilter.

Definition at line 19 of file PFTauVertexSelector.cc.

References beamSpotSrc_, reco::TrackBase::dz(), dZ_, filterOnNTaus_, i, edm::HandleBase::isValid(), recoCandidateSrc_, tauSrc_, testEve_cfg::tracks, trackSrc_, trigger::TriggerElectron, triggerFilterElectronsSrc_, triggerFilterMuonsSrc_, trigger::TriggerMuon, useBeamSpot_, useLeadingRecoCandidate_, useLeadingTrack_, useTriggerFilterElectrons_, useTriggerFilterMuons_, useVertex_, and vertexSrc_.

                                                                                 {

   math::XYZPoint vertexPoint;
   bool vertexAvailable=false;
   
   if(useBeamSpot_)
   {
       edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
       event.getByLabel(beamSpotSrc_,recoBeamSpotHandle);
       if (recoBeamSpotHandle.isValid()){
         vertexPoint = recoBeamSpotHandle->position();
         vertexAvailable = true;
       }
   }
 
   if(useVertex_)
   {
       edm::Handle<edm::View<reco::Vertex> > recoVertexHandle;
       event.getByLabel(vertexSrc_,recoVertexHandle);
       if ((recoVertexHandle.isValid()) && (recoVertexHandle->size()>0)){
         vertexPoint = recoVertexHandle->at(0).position();
         vertexAvailable = true;
       }
   }
 
   const reco::Track* track=0;
   double maxpt=0.;
   
   if (useLeadingTrack_)
   {
       edm::Handle<edm::View<reco::Track> > tracks;
       for( std::vector<edm::InputTag>::const_iterator trackSrc = trackSrc_.begin(); trackSrc != trackSrc_.end(); ++ trackSrc ) {
           event.getByLabel(*trackSrc, tracks);
           if ((tracks.isValid())&&(tracks->size()>0)){
               for (unsigned i = 0; i < tracks->size(); ++i) {
                  double pt=tracks->ptrAt(i)->pt();
                  if(pt>maxpt)
                  {
                      track = &*tracks->ptrAt(i);
                      maxpt=pt;
                  }
               }
           }
       }
   }
   
   if (useLeadingRecoCandidate_)
   {
       edm::Handle<edm::View<reco::RecoCandidate> > recocandidates;
       for( std::vector<edm::InputTag>::const_iterator recoCandidateSrc = recoCandidateSrc_.begin(); recoCandidateSrc != recoCandidateSrc_.end(); ++ recoCandidateSrc ) {
           event.getByLabel(*recoCandidateSrc, recocandidates);
           if ((recocandidates.isValid())&&(recocandidates->size()>0)){
               for (unsigned i = 0; i < recocandidates->size(); ++i) {
                  double pt=recocandidates->ptrAt(i)->pt();
                  if(pt>maxpt)
                  {
                      track = dynamic_cast<const reco::Track*>(recocandidates->ptrAt(i)->bestTrack());
                      maxpt=pt;
                  }
               }
           }
       }
   }
   
   if (useTriggerFilterElectrons_)
   {
       edm::Handle<trigger::TriggerFilterObjectWithRefs> triggerfilter;
       event.getByLabel(triggerFilterElectronsSrc_, triggerfilter);
       std::vector<reco::ElectronRef> recocandidates;
       triggerfilter->getObjects(trigger::TriggerElectron,recocandidates);
       if ((recocandidates.size()>0)){
           for (unsigned i = 0; i < recocandidates.size(); ++i) {
              double pt=recocandidates.at(i)->pt();
              if(pt>maxpt)
              {
                  track = dynamic_cast<const reco::Track*>(recocandidates.at(i)->bestTrack());
                  maxpt=pt;
              }
           }
       }
   }
   
   if (useTriggerFilterMuons_)
   {
       edm::Handle<trigger::TriggerFilterObjectWithRefs> triggerfilter;
       event.getByLabel(triggerFilterMuonsSrc_, triggerfilter);
       std::vector<reco::RecoChargedCandidateRef> recocandidates;
       triggerfilter->getObjects(trigger::TriggerMuon,recocandidates);
       if ((recocandidates.size()>0)){
           for (unsigned i = 0; i < recocandidates.size(); ++i) {
              double pt=recocandidates.at(i)->pt();
              if(pt>maxpt)
              {
                  track = dynamic_cast<const reco::Track*>(recocandidates.at(i)->bestTrack());
                  maxpt=pt;
              }
           }
       }
   }
   
   reco::PFTauCollection* selTaus = new reco::PFTauCollection;
   edm::Handle<edm::View<reco::PFTau> > taus;
   event.getByLabel(tauSrc_, taus);
   for( edm::View<reco::PFTau>::const_iterator pfTau = taus->begin(); pfTau != taus->end(); ++ pfTau ) {
       // if no leading track assigned skip
       if ((!pfTau->leadPFChargedHadrCand().isNonnull())||
           (!pfTau->leadPFChargedHadrCand()->trackRef().isNonnull()))
          continue;

       if(vertexAvailable)
       {
           // select by z position of leading track at vertex
           if ((useLeadingTrack_)||(useLeadingRecoCandidate_)||(useTriggerFilterElectrons_)||(useTriggerFilterMuons_))
           {
               if((track)&&(fabs(pfTau->leadPFChargedHadrCand()->trackRef()->dz(vertexPoint) - track->dz(vertexPoint)) < dZ_))
                   selTaus->push_back(*pfTau);
           }
           // select by z position of leading vertex
           else
           {
               if (fabs(pfTau->leadPFChargedHadrCand()->trackRef()->dz(vertexPoint))<dZ_)
                   selTaus->push_back(*pfTau);
           }
       }
       else
       {
           // select by z position of leading track at (0,0,0)
           if ((useLeadingTrack_)||(useLeadingRecoCandidate_)||(useTriggerFilterElectrons_)||(useTriggerFilterMuons_))
           {
               if((track)&&(fabs(pfTau->leadPFChargedHadrCand()->trackRef()->dz() - track->dz()) < dZ_))
                   selTaus->push_back(*pfTau);
           }
       }
   }
   unsigned filterTaus=selTaus->size();
   std::auto_ptr<reco::PFTauCollection> selectedTaus(selTaus);
   event.put(selectedTaus);
   
   return (filterTaus>=filterOnNTaus_);
}

Member Data Documentation

Definition at line 50 of file PFTauVertexSelector.h.

Referenced by filter(), and PFTauVertexSelector().

double PFTauVertexSelector::dZ_ [private]

Definition at line 59 of file PFTauVertexSelector.h.

Referenced by filter(), and PFTauVertexSelector().

Definition at line 60 of file PFTauVertexSelector.h.

Referenced by filter(), and PFTauVertexSelector().

Definition at line 54 of file PFTauVertexSelector.h.

Referenced by filter(), and PFTauVertexSelector().

Definition at line 46 of file PFTauVertexSelector.h.

Referenced by filter(), and PFTauVertexSelector().

Definition at line 52 of file PFTauVertexSelector.h.

Referenced by filter(), and PFTauVertexSelector().

Definition at line 56 of file PFTauVertexSelector.h.

Referenced by filter(), and PFTauVertexSelector().

Definition at line 58 of file PFTauVertexSelector.h.

Referenced by filter(), and PFTauVertexSelector().

Definition at line 49 of file PFTauVertexSelector.h.

Referenced by filter(), and PFTauVertexSelector().

Definition at line 53 of file PFTauVertexSelector.h.

Referenced by filter(), and PFTauVertexSelector().

Definition at line 51 of file PFTauVertexSelector.h.

Referenced by filter(), and PFTauVertexSelector().

Definition at line 55 of file PFTauVertexSelector.h.

Referenced by filter(), and PFTauVertexSelector().

Definition at line 57 of file PFTauVertexSelector.h.

Referenced by filter(), and PFTauVertexSelector().

Definition at line 47 of file PFTauVertexSelector.h.

Referenced by filter(), and PFTauVertexSelector().

Definition at line 48 of file PFTauVertexSelector.h.

Referenced by filter(), and PFTauVertexSelector().