CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

PFRecoTauDiscriminationByFlightPathSignificance Class Reference

Inheritance diagram for PFRecoTauDiscriminationByFlightPathSignificance:
TauDiscriminationProducerBase< TauType, TauDiscriminator > edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

List of all members.

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 TransientTrackBuildertransientTrackBuilder
reco::tau::RecoTauVertexAssociator vertexAssociator_
bool withPVError

Detailed Description

Definition at line 25 of file PFRecoTauDiscriminationByFlightPathSignificance.cc.


Constructor & Destructor Documentation

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]

Member Function Documentation

void PFRecoTauDiscriminationByFlightPathSignificance::beginEvent ( const edm::Event iEvent,
const edm::EventSetup iSetup 
) [virtual]
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();
}

Member Data Documentation