Public Member Functions | |
void | beginEvent (const edm::Event &, const edm::EventSetup &) |
double | discriminate (const reco::PFTauRef &) |
PFRecoTauDiscriminationByFlightPathSignificance (const ParameterSet &iConfig) | |
~PFRecoTauDiscriminationByFlightPathSignificance () | |
Private Member Functions | |
double | threeProngFlightPathSig (const PFTauRef &) |
double | vertexSignificance (reco::Vertex &, reco::Vertex &, GlobalVector &) |
Private Attributes | |
bool | booleanOutput |
double | flightPathSig |
const TransientTrackBuilder * | transientTrackBuilder |
reco::tau::RecoTauVertexAssociator | vertexAssociator_ |
bool | withPVError |
Definition at line 25 of file PFRecoTauDiscriminationByFlightPathSignificance.cc.
PFRecoTauDiscriminationByFlightPathSignificance::PFRecoTauDiscriminationByFlightPathSignificance | ( | const ParameterSet & | iConfig | ) | [inline, explicit] |
Definition at line 28 of file PFRecoTauDiscriminationByFlightPathSignificance.cc.
References edm::ParameterSet::getParameter().
:PFTauDiscriminationProducerBase(iConfig), vertexAssociator_(iConfig.getParameter<ParameterSet>("qualityCuts")) { flightPathSig = iConfig.getParameter<double>("flightPathSig"); withPVError = iConfig.getParameter<bool>("UsePVerror"); booleanOutput = iConfig.getParameter<bool>("BooleanOutput"); }
PFRecoTauDiscriminationByFlightPathSignificance::~PFRecoTauDiscriminationByFlightPathSignificance | ( | ) | [inline] |
Definition at line 37 of file PFRecoTauDiscriminationByFlightPathSignificance.cc.
{}
void PFRecoTauDiscriminationByFlightPathSignificance::beginEvent | ( | const edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
Reimplemented from TauDiscriminationProducerBase< TauType, TauDiscriminator >.
Definition at line 55 of file PFRecoTauDiscriminationByFlightPathSignificance.cc.
References edm::EventSetup::get(), and edm::ESHandle< T >::product().
{ vertexAssociator_.setEvent(iEvent); // Transient Tracks edm::ESHandle<TransientTrackBuilder> builder; iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",builder); transientTrackBuilder = builder.product(); }
double PFRecoTauDiscriminationByFlightPathSignificance::discriminate | ( | const reco::PFTauRef & | tau | ) |
Definition at line 67 of file PFRecoTauDiscriminationByFlightPathSignificance.cc.
{ if(booleanOutput) return ( threeProngFlightPathSig(tau) > flightPathSig ? 1. : 0. ); return threeProngFlightPathSig(tau); }
double PFRecoTauDiscriminationByFlightPathSignificance::threeProngFlightPathSig | ( | const PFTauRef & | tau | ) | [private] |
Definition at line 73 of file PFRecoTauDiscriminationByFlightPathSignificance.cc.
References edm::RefVector< C, T, F >::begin(), edm::RefVector< C, T, F >::end(), edm::Ref< C, T, F >::isNonnull(), edm::Ref< C, T, F >::isNull(), TransientVertex::isValid(), reco::PFCandidate::trackRef(), and KalmanVertexFitter::vertex().
{ double flightPathSignificance = 0; reco::VertexRef primaryVertex = vertexAssociator_.associatedVertex(*tau); if (primaryVertex.isNull()) { edm::LogError("FlightPathSignficance") << "Could not get vertex associated" << " to tau, returning -999!" << std::endl; return -999; } //Secondary vertex const PFCandidateRefVector pfSignalCandidates = tau->signalPFChargedHadrCands(); vector<TransientTrack> transientTracks; RefVector<PFCandidateCollection>::const_iterator iTrack; for(iTrack = pfSignalCandidates.begin(); iTrack!= pfSignalCandidates.end(); iTrack++){ const PFCandidate& pfCand = *(iTrack->get()); if(pfCand.trackRef().isNonnull()){ const TransientTrack transientTrack = transientTrackBuilder->build(pfCand.trackRef()); transientTracks.push_back(transientTrack); } } if(transientTracks.size() > 1){ KalmanVertexFitter kvf(true); TransientVertex tv = kvf.vertex(transientTracks); if(tv.isValid()){ GlobalVector tauDir(tau->px(), tau->py(), tau->pz()); Vertex secVer = tv; // We have to un-const the PV for some reason reco::Vertex primaryVertexNonConst = *primaryVertex; flightPathSignificance = vertexSignificance(primaryVertexNonConst,secVer,tauDir); } } return flightPathSignificance; }
double PFRecoTauDiscriminationByFlightPathSignificance::vertexSignificance | ( | reco::Vertex & | pv, |
reco::Vertex & | sv, | ||
GlobalVector & | direction | ||
) | [private] |
Definition at line 113 of file PFRecoTauDiscriminationByFlightPathSignificance.cc.
References reco::SecondaryVertex::computeDist3d(), and Measurement1D::significance().
{ return SecondaryVertex::computeDist3d(pv,sv,direction,withPVError).significance(); }
Definition at line 48 of file PFRecoTauDiscriminationByFlightPathSignificance.cc.
double PFRecoTauDiscriminationByFlightPathSignificance::flightPathSig [private] |
Definition at line 49 of file PFRecoTauDiscriminationByFlightPathSignificance.cc.
const TransientTrackBuilder* PFRecoTauDiscriminationByFlightPathSignificance::transientTrackBuilder [private] |
Definition at line 52 of file PFRecoTauDiscriminationByFlightPathSignificance.cc.
reco::tau::RecoTauVertexAssociator PFRecoTauDiscriminationByFlightPathSignificance::vertexAssociator_ [private] |
Definition at line 46 of file PFRecoTauDiscriminationByFlightPathSignificance.cc.
bool PFRecoTauDiscriminationByFlightPathSignificance::withPVError [private] |
Definition at line 50 of file PFRecoTauDiscriminationByFlightPathSignificance.cc.