CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/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 : public PFTauDiscriminationProducerBase  {
00025     public:
00026         explicit PFRecoTauDiscriminationByFlightPathSignificance(const ParameterSet& iConfig):PFTauDiscriminationProducerBase(iConfig),
00027                                                                                qualityCuts_(iConfig.getParameter<ParameterSet>("qualityCuts")){  // retrieve quality cuts
00028                 flightPathSig           = iConfig.getParameter<double>("flightPathSig");
00029                 withPVError             = iConfig.getParameter<bool>("UsePVerror");
00030 
00031                 PVProducer              = iConfig.getParameter<edm::InputTag>("PVProducer");
00032 
00033                 booleanOutput           = iConfig.getParameter<bool>("BooleanOutput");
00034         }
00035 
00036         ~PFRecoTauDiscriminationByFlightPathSignificance(){}
00037 
00038         void beginEvent(const edm::Event&, const edm::EventSetup&);
00039         double discriminate(const reco::PFTauRef&);
00040 
00041     private:
00042         double threeProngFlightPathSig(const PFTauRef&);
00043         double vertexSignificance(reco::Vertex&,reco::Vertex&,GlobalVector&);
00044 
00045         PFTauQualityCutWrapper qualityCuts_;
00046 
00047         double flightPathSig;
00048         bool withPVError;
00049 
00050         reco::Vertex primaryVertex;
00051         const TransientTrackBuilder* transientTrackBuilder;
00052         edm::InputTag PVProducer;
00053 
00054         bool booleanOutput;
00055 };
00056 
00057 void PFRecoTauDiscriminationByFlightPathSignificance::beginEvent(const Event& iEvent, const EventSetup& iSetup){
00058 
00059 //Primary vertex
00060         edm::Handle<edm::View<reco::Vertex> > vertexHandle;
00061         iEvent.getByLabel(PVProducer, vertexHandle);
00062         const edm::View<reco::Vertex>& vertexCollection(*vertexHandle);
00063 
00064         primaryVertex = *(vertexCollection.begin());
00065 
00066         // Transient Tracks
00067         edm::ESHandle<TransientTrackBuilder> builder;
00068         iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",builder);
00069         transientTrackBuilder = builder.product();
00070 
00071 }
00072 
00073 double PFRecoTauDiscriminationByFlightPathSignificance::discriminate(const PFTauRef& tau){
00074 
00075         if(booleanOutput) return ( threeProngFlightPathSig(tau) > flightPathSig ? 1. : 0. );
00076         return threeProngFlightPathSig(tau);
00077 }
00078 
00079 double PFRecoTauDiscriminationByFlightPathSignificance::threeProngFlightPathSig(const PFTauRef& tau){
00080         double flightPathSignificance = 0;
00081 
00082 //Secondary vertex
00083         const PFCandidateRefVector pfSignalCandidates = tau->signalPFChargedHadrCands();
00084         vector<TransientTrack> transientTracks;
00085         RefVector<PFCandidateCollection>::const_iterator iTrack;
00086         for(iTrack = pfSignalCandidates.begin(); iTrack!= pfSignalCandidates.end(); iTrack++){
00087                 const PFCandidate& pfCand = *(iTrack->get());
00088                 if(pfCand.trackRef().isNonnull()){
00089                   const TransientTrack transientTrack = transientTrackBuilder->build(pfCand.trackRef());
00090                   transientTracks.push_back(transientTrack);
00091                 }
00092         }
00093         if(transientTracks.size() > 1){
00094                 KalmanVertexFitter kvf(true);
00095                 TransientVertex tv = kvf.vertex(transientTracks);
00096 
00097                 if(tv.isValid()){
00098                         GlobalVector tauDir(tau->px(),
00099                                             tau->py(),
00100                                             tau->pz());
00101                         Vertex secVer = tv;
00102                         flightPathSignificance = vertexSignificance(primaryVertex,secVer,tauDir);
00103                 }
00104         }
00105         return flightPathSignificance;
00106 }
00107 
00108 double PFRecoTauDiscriminationByFlightPathSignificance::vertexSignificance(reco::Vertex& pv, Vertex& sv,GlobalVector& direction){
00109         return SecondaryVertex::computeDist3d(pv,sv,direction,withPVError).significance();
00110 }
00111 
00112 DEFINE_FWK_MODULE(PFRecoTauDiscriminationByFlightPathSignificance);
00113