CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoTauTag/RecoTau/plugins/PFRecoTauDiscriminationByFlightPathSignificance.cc

Go to the documentation of this file.
00001 #include "RecoTauTag/RecoTau/interface/TauDiscriminationProducerBase.h"
00002 #include "RecoTauTag/TauTagTools/interface/PFTauQualityCutWrapper.h"
00003 #include "FWCore/Utilities/interface/InputTag.h"
00004 
00005 /* class PFRecoTauDiscriminationByFlightPathSignificance
00006  * created : August 30 2010,
00007  * contributors : Sami Lehti (sami.lehti@cern.ch ; HIP, Helsinki)
00008  * based on H+ tau ID by Lauri Wendland
00009  */
00010 
00011 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00012 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
00013 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00014 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexFitter.h"
00015 #include "RecoVertex/VertexPrimitives/interface/TransientVertex.h"
00016 #include "RecoBTag/SecondaryVertex/interface/SecondaryVertex.h"
00017 
00018 #include "TLorentzVector.h"
00019 
00020 using namespace reco;
00021 using namespace std;
00022 using namespace edm;
00023 
00024 class PFRecoTauDiscriminationByFlightPathSignificance
00025 : public PFTauDiscriminationProducerBase  {
00026   public:
00027     explicit PFRecoTauDiscriminationByFlightPathSignificance(const ParameterSet& iConfig):PFTauDiscriminationProducerBase(iConfig),
00028     qualityCuts_(iConfig.getParameter<ParameterSet>("qualityCuts")){  // retrieve quality cuts
00029       flightPathSig             = iConfig.getParameter<double>("flightPathSig");
00030       withPVError               = iConfig.getParameter<bool>("UsePVerror");
00031 
00032       PVProducer                = iConfig.getParameter<edm::InputTag>("PVProducer");
00033 
00034       booleanOutput             = iConfig.getParameter<bool>("BooleanOutput");
00035     }
00036 
00037     ~PFRecoTauDiscriminationByFlightPathSignificance(){}
00038 
00039     void beginEvent(const edm::Event&, const edm::EventSetup&);
00040     double discriminate(const reco::PFTauRef&);
00041 
00042   private:
00043     double threeProngFlightPathSig(const PFTauRef&);
00044     double vertexSignificance(reco::Vertex&,reco::Vertex&,GlobalVector&);
00045 
00046     PFTauQualityCutWrapper qualityCuts_;
00047 
00048     double flightPathSig;
00049     bool withPVError;
00050 
00051     reco::Vertex primaryVertex;
00052     const TransientTrackBuilder* transientTrackBuilder;
00053     edm::InputTag PVProducer;
00054 
00055     bool booleanOutput;
00056 };
00057 
00058 void PFRecoTauDiscriminationByFlightPathSignificance::beginEvent(const Event& iEvent, const EventSetup& iSetup){
00059 
00060   //Primary vertex
00061   edm::Handle<edm::View<reco::Vertex> > vertexHandle;
00062   iEvent.getByLabel(PVProducer, vertexHandle);
00063   const edm::View<reco::Vertex>& vertexCollection(*vertexHandle);
00064 
00065   primaryVertex = *(vertexCollection.begin());
00066 
00067   // Transient Tracks
00068   edm::ESHandle<TransientTrackBuilder> builder;
00069   iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",builder);
00070   transientTrackBuilder = builder.product();
00071 
00072 }
00073 
00074 double PFRecoTauDiscriminationByFlightPathSignificance::discriminate(const PFTauRef& tau){
00075 
00076   if(booleanOutput) return ( threeProngFlightPathSig(tau) > flightPathSig ? 1. : 0. );
00077   return threeProngFlightPathSig(tau);
00078 }
00079 
00080 double PFRecoTauDiscriminationByFlightPathSignificance::threeProngFlightPathSig(
00081     const PFTauRef& tau){
00082   double flightPathSignificance = 0;
00083 
00084   //Secondary vertex
00085   const PFCandidateRefVector pfSignalCandidates = tau->signalPFChargedHadrCands();
00086   vector<TransientTrack> transientTracks;
00087   RefVector<PFCandidateCollection>::const_iterator iTrack;
00088   for(iTrack = pfSignalCandidates.begin(); iTrack!= pfSignalCandidates.end(); iTrack++){
00089     const PFCandidate& pfCand = *(iTrack->get());
00090     if(pfCand.trackRef().isNonnull()){
00091       const TransientTrack transientTrack = transientTrackBuilder->build(pfCand.trackRef());
00092       transientTracks.push_back(transientTrack);
00093     }
00094   }
00095   if(transientTracks.size() > 1){
00096     KalmanVertexFitter kvf(true);
00097     TransientVertex tv = kvf.vertex(transientTracks);
00098 
00099     if(tv.isValid()){
00100       GlobalVector tauDir(tau->px(),
00101           tau->py(),
00102           tau->pz());
00103       Vertex secVer = tv;
00104       flightPathSignificance = vertexSignificance(primaryVertex,secVer,tauDir);
00105     }
00106   }
00107   return flightPathSignificance;
00108 }
00109 
00110 double PFRecoTauDiscriminationByFlightPathSignificance::vertexSignificance(
00111     reco::Vertex& pv, Vertex& sv,GlobalVector& direction){
00112   return SecondaryVertex::computeDist3d(pv,sv,direction,withPVError).significance();
00113 }
00114 
00115 DEFINE_FWK_MODULE(PFRecoTauDiscriminationByFlightPathSignificance);
00116