Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Updated PFRecoTauDiscriminationByNProngs #3540

Merged
merged 1 commit into from Apr 30, 2014
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
55 changes: 43 additions & 12 deletions RecoTauTag/RecoTau/plugins/PFRecoTauDiscriminationByNProngs.cc
@@ -1,10 +1,15 @@
#include "RecoTauTag/RecoTau/interface/TauDiscriminationProducerBase.h"
#include "FWCore/Utilities/interface/InputTag.h"
#include "RecoTauTag/RecoTau/interface/RecoTauQualityCuts.h"
#include "RecoTauTag/RecoTau/interface/RecoTauVertexAssociator.h"
#include <boost/foreach.hpp>
#include "DataFormats/VertexReco/interface/Vertex.h"

/* class PFRecoTauDiscriminationByNProngs
* created : August 30 2010,
* contributors : Sami Lehti (sami.lehti@cern.ch ; HIP, Helsinki)
* based on H+ tau ID by Lauri Wendland
* Modified April 16 2014 by S.Lehti
*/

using namespace reco;
Expand All @@ -13,36 +18,62 @@ using namespace edm;

class PFRecoTauDiscriminationByNProngs : public PFTauDiscriminationProducerBase {
public:
explicit PFRecoTauDiscriminationByNProngs(const ParameterSet& iConfig):PFTauDiscriminationProducerBase(iConfig){
nprongs = iConfig.getParameter<uint32_t>("nProngs");
booleanOutput = iConfig.getParameter<bool>("BooleanOutput");
}

explicit PFRecoTauDiscriminationByNProngs(const ParameterSet&);
~PFRecoTauDiscriminationByNProngs(){}

void beginEvent(const edm::Event&, const edm::EventSetup&) override;
double discriminate(const reco::PFTauRef&) override;

private:
std::auto_ptr<tau::RecoTauQualityCuts> qcuts_;
std::auto_ptr<tau::RecoTauVertexAssociator> vertexAssociator_;

uint32_t nprongs;
uint32_t minN,maxN;
bool booleanOutput;
edm::ParameterSet qualityCuts;
};

void PFRecoTauDiscriminationByNProngs::beginEvent(const Event& iEvent, const EventSetup& iSetup){}
PFRecoTauDiscriminationByNProngs::PFRecoTauDiscriminationByNProngs(const ParameterSet& iConfig):
PFTauDiscriminationProducerBase(iConfig),
qualityCuts(iConfig.getParameterSet("qualityCuts"))
{
minN = iConfig.getParameter<uint32_t>("MinN");
maxN = iConfig.getParameter<uint32_t>("MaxN");
booleanOutput = iConfig.getParameter<bool>("BooleanOutput");

qcuts_.reset(new tau::RecoTauQualityCuts(qualityCuts.getParameterSet("signalQualityCuts")));
vertexAssociator_.reset(new tau::RecoTauVertexAssociator(qualityCuts,consumesCollector()));
}

void PFRecoTauDiscriminationByNProngs::beginEvent(const Event& iEvent, const EventSetup& iSetup){
vertexAssociator_->setEvent(iEvent);
}

double PFRecoTauDiscriminationByNProngs::discriminate(const PFTauRef& tau){

bool accepted = false;
int np = tau->signalTracks().size();
reco::VertexRef pv = vertexAssociator_->associatedVertex(*tau);
const PFCandidatePtr leadingTrack = tau->leadPFChargedHadrCand();

if((np == 1 && (nprongs == 1 || nprongs == 0)) ||
(np == 3 && (nprongs == 3 || nprongs == 0)) ) accepted = true;
uint np = 0;
if(leadingTrack.isNonnull() && pv.isNonnull()){
qcuts_->setPV(pv);
qcuts_->setLeadTrack(tau->leadPFChargedHadrCand());

BOOST_FOREACH( const reco::PFCandidatePtr& cand, tau->signalPFChargedHadrCands() ) {
if ( qcuts_->filterCandRef(cand) ) np++;
}
}

bool accepted = false;
if(maxN == 0){
if(np == 1 || np == 3) accepted = true;
}else{
if(np >= minN && np <= maxN) accepted = true;
}

if(!accepted) np = 0;
if(booleanOutput) return accepted;
return np;
}

DEFINE_FWK_MODULE(PFRecoTauDiscriminationByNProngs);

14 changes: 10 additions & 4 deletions RecoTauTag/RecoTau/python/PFRecoTauDiscriminationByNProngs_cfi.py
@@ -1,13 +1,19 @@
import FWCore.ParameterSet.Config as cms

from RecoTauTag.RecoTau.PFRecoTauQualityCuts_cfi import *
from RecoTauTag.RecoTau.TauDiscriminatorTools import requireLeadTrack
#from RecoTauTag.RecoTau.TauDiscriminatorTools import requireLeadTrack
from RecoTauTag.RecoTau.TauDiscriminatorTools import noPrediscriminants
from RecoTauTag.RecoTau.PFRecoTauQualityCuts_cfi import PFTauQualityCuts

pfRecoTauDiscriminationByNProngs = cms.EDProducer("PFRecoTauDiscriminationByNProngs",
PFTauProducer = cms.InputTag('pfRecoTauProducer'), #tau collection to discriminate
PFTauProducer = cms.InputTag('combinatoricRecoTaus'),

Prediscriminants = requireLeadTrack,
# Prediscriminants = requireLeadTrack,
Prediscriminants = noPrediscriminants,
BooleanOutput = cms.bool(True),

nProngs = cms.uint32(0), # number of prongs required: 0=1||3, 1, 3
MaxN = cms.uint32(0), # number of prongs required: 0=1||3
MinN = cms.uint32(1),

qualityCuts = PFTauQualityCuts
)