CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/RecoTauTag/HLTProducers/src/PFTauVertexSelector.cc

Go to the documentation of this file.
00001 #include "RecoTauTag/HLTProducers/interface/PFTauVertexSelector.h"
00002 #include "DataFormats/VertexReco/interface/Vertex.h"
00003 #include "DataFormats/TrackReco/interface/Track.h"
00004 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
00005 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
00006 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
00007 #include "DataFormats/EgammaCandidates/interface/Electron.h"
00008 #include "DataFormats/Math/interface/Point3D.h"
00009 #include "DataFormats/Math/interface/Error.h"
00010 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
00011 
00012 /* 
00013  * class PFTauVertexSelector
00014  * created : January 26 2012,
00015  * revised : Wed Jan 26 11:13:04 PDT 2012
00016  * Authors : Andreas Hinzmann (CERN)
00017  */
00018 
00019 bool PFTauVertexSelector::filter(edm::Event& event, const edm::EventSetup& eventSetup) {
00020 
00021    math::XYZPoint vertexPoint;
00022    bool vertexAvailable=false;
00023    
00024    if(useBeamSpot_)
00025    {
00026        edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00027        event.getByLabel(beamSpotSrc_,recoBeamSpotHandle);
00028        if (recoBeamSpotHandle.isValid()){
00029          vertexPoint = recoBeamSpotHandle->position();
00030          vertexAvailable = true;
00031        }
00032    }
00033  
00034    if(useVertex_)
00035    {
00036        edm::Handle<edm::View<reco::Vertex> > recoVertexHandle;
00037        event.getByLabel(vertexSrc_,recoVertexHandle);
00038        if ((recoVertexHandle.isValid()) && (recoVertexHandle->size()>0)){
00039          vertexPoint = recoVertexHandle->at(0).position();
00040          vertexAvailable = true;
00041        }
00042    }
00043  
00044    const reco::Track* track=0;
00045    double maxpt=0.;
00046    
00047    if (useLeadingTrack_)
00048    {
00049        edm::Handle<edm::View<reco::Track> > tracks;
00050        for( std::vector<edm::InputTag>::const_iterator trackSrc = trackSrc_.begin(); trackSrc != trackSrc_.end(); ++ trackSrc ) {
00051            event.getByLabel(*trackSrc, tracks);
00052            if ((tracks.isValid())&&(tracks->size()>0)){
00053                for (unsigned i = 0; i < tracks->size(); ++i) {
00054                   double pt=tracks->ptrAt(i)->pt();
00055                   if(pt>maxpt)
00056                   {
00057                       track = &*tracks->ptrAt(i);
00058                       maxpt=pt;
00059                   }
00060                }
00061            }
00062        }
00063    }
00064    
00065    if (useLeadingRecoCandidate_)
00066    {
00067        edm::Handle<edm::View<reco::RecoCandidate> > recocandidates;
00068        for( std::vector<edm::InputTag>::const_iterator recoCandidateSrc = recoCandidateSrc_.begin(); recoCandidateSrc != recoCandidateSrc_.end(); ++ recoCandidateSrc ) {
00069            event.getByLabel(*recoCandidateSrc, recocandidates);
00070            if ((recocandidates.isValid())&&(recocandidates->size()>0)){
00071                for (unsigned i = 0; i < recocandidates->size(); ++i) {
00072                   double pt=recocandidates->ptrAt(i)->pt();
00073                   if(pt>maxpt)
00074                   {
00075                       track = dynamic_cast<const reco::Track*>(recocandidates->ptrAt(i)->bestTrack());
00076                       maxpt=pt;
00077                   }
00078                }
00079            }
00080        }
00081    }
00082    
00083    if (useTriggerFilterElectrons_)
00084    {
00085        edm::Handle<trigger::TriggerFilterObjectWithRefs> triggerfilter;
00086        event.getByLabel(triggerFilterElectronsSrc_, triggerfilter);
00087        std::vector<reco::ElectronRef> recocandidates;
00088        triggerfilter->getObjects(trigger::TriggerElectron,recocandidates);
00089        if ((recocandidates.size()>0)){
00090            for (unsigned i = 0; i < recocandidates.size(); ++i) {
00091               double pt=recocandidates.at(i)->pt();
00092               if(pt>maxpt)
00093               {
00094                   track = dynamic_cast<const reco::Track*>(recocandidates.at(i)->bestTrack());
00095                   maxpt=pt;
00096               }
00097            }
00098        }
00099    }
00100    
00101    if (useTriggerFilterMuons_)
00102    {
00103        edm::Handle<trigger::TriggerFilterObjectWithRefs> triggerfilter;
00104        event.getByLabel(triggerFilterMuonsSrc_, triggerfilter);
00105        std::vector<reco::RecoChargedCandidateRef> recocandidates;
00106        triggerfilter->getObjects(trigger::TriggerMuon,recocandidates);
00107        if ((recocandidates.size()>0)){
00108            for (unsigned i = 0; i < recocandidates.size(); ++i) {
00109               double pt=recocandidates.at(i)->pt();
00110               if(pt>maxpt)
00111               {
00112                   track = dynamic_cast<const reco::Track*>(recocandidates.at(i)->bestTrack());
00113                   maxpt=pt;
00114               }
00115            }
00116        }
00117    }
00118    
00119    reco::PFTauCollection* selTaus = new reco::PFTauCollection;
00120    edm::Handle<edm::View<reco::PFTau> > taus;
00121    event.getByLabel(tauSrc_, taus);
00122    for( edm::View<reco::PFTau>::const_iterator pfTau = taus->begin(); pfTau != taus->end(); ++ pfTau ) {
00123        // if no leading track assigned skip
00124        if ((!pfTau->leadPFChargedHadrCand().isNonnull())||
00125            (!pfTau->leadPFChargedHadrCand()->trackRef().isNonnull()))
00126           continue;
00127 
00128        if(vertexAvailable)
00129        {
00130            // select by z position of leading track at vertex
00131            if ((useLeadingTrack_)||(useLeadingRecoCandidate_)||(useTriggerFilterElectrons_)||(useTriggerFilterMuons_))
00132            {
00133                if((track)&&(fabs(pfTau->leadPFChargedHadrCand()->trackRef()->dz(vertexPoint) - track->dz(vertexPoint)) < dZ_))
00134                    selTaus->push_back(*pfTau);
00135            }
00136            // select by z position of leading vertex
00137            else
00138            {
00139                if (fabs(pfTau->leadPFChargedHadrCand()->trackRef()->dz(vertexPoint))<dZ_)
00140                    selTaus->push_back(*pfTau);
00141            }
00142        }
00143        else
00144        {
00145            // select by z position of leading track at (0,0,0)
00146            if ((useLeadingTrack_)||(useLeadingRecoCandidate_)||(useTriggerFilterElectrons_)||(useTriggerFilterMuons_))
00147            {
00148                if((track)&&(fabs(pfTau->leadPFChargedHadrCand()->trackRef()->dz() - track->dz()) < dZ_))
00149                    selTaus->push_back(*pfTau);
00150            }
00151        }
00152    }
00153    unsigned filterTaus=selTaus->size();
00154    std::auto_ptr<reco::PFTauCollection> selectedTaus(selTaus);
00155    event.put(selectedTaus);
00156    
00157    return (filterTaus>=filterOnNTaus_);
00158 }