CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_1/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 #include "RecoTauTag/RecoTau/interface/RecoTauVertexAssociator.h"
00018 
00019 #include "TLorentzVector.h"
00020 
00021 using namespace reco;
00022 using namespace std;
00023 using namespace edm;
00024 
00025 class PFRecoTauDiscriminationByFlightPathSignificance
00026   : public PFTauDiscriminationProducerBase  {
00027   public:
00028     explicit PFRecoTauDiscriminationByFlightPathSignificance(
00029         const ParameterSet& iConfig)
00030       :PFTauDiscriminationProducerBase(iConfig),
00031       vertexAssociator_(iConfig.getParameter<ParameterSet>("qualityCuts")) {
00032       flightPathSig             = iConfig.getParameter<double>("flightPathSig");
00033       withPVError               = iConfig.getParameter<bool>("UsePVerror");
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     reco::tau::RecoTauVertexAssociator vertexAssociator_;
00047 
00048     bool booleanOutput;
00049     double flightPathSig;
00050     bool withPVError;
00051 
00052     const TransientTrackBuilder* transientTrackBuilder;
00053 };
00054 
00055 void PFRecoTauDiscriminationByFlightPathSignificance::beginEvent(
00056     const Event& iEvent, const EventSetup& iSetup){
00057 
00058   vertexAssociator_.setEvent(iEvent);
00059 
00060   // Transient Tracks
00061   edm::ESHandle<TransientTrackBuilder> builder;
00062   iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",builder);
00063   transientTrackBuilder = builder.product();
00064 
00065 }
00066 
00067 double PFRecoTauDiscriminationByFlightPathSignificance::discriminate(const PFTauRef& tau){
00068 
00069   if(booleanOutput) return ( threeProngFlightPathSig(tau) > flightPathSig ? 1. : 0. );
00070   return threeProngFlightPathSig(tau);
00071 }
00072 
00073 double PFRecoTauDiscriminationByFlightPathSignificance::threeProngFlightPathSig(
00074     const PFTauRef& tau){
00075   double flightPathSignificance = 0;
00076 
00077   reco::VertexRef primaryVertex = vertexAssociator_.associatedVertex(*tau);
00078 
00079   if (primaryVertex.isNull()) {
00080     edm::LogError("FlightPathSignficance") << "Could not get vertex associated"
00081       << " to tau, returning -999!" << std::endl;
00082     return -999;
00083   }
00084 
00085   //Secondary vertex
00086   const PFCandidateRefVector pfSignalCandidates = tau->signalPFChargedHadrCands();
00087   vector<TransientTrack> transientTracks;
00088   RefVector<PFCandidateCollection>::const_iterator iTrack;
00089   for(iTrack = pfSignalCandidates.begin(); iTrack!= pfSignalCandidates.end(); iTrack++){
00090     const PFCandidate& pfCand = *(iTrack->get());
00091     if(pfCand.trackRef().isNonnull()){
00092       const TransientTrack transientTrack = transientTrackBuilder->build(pfCand.trackRef());
00093       transientTracks.push_back(transientTrack);
00094     }
00095   }
00096   if(transientTracks.size() > 1){
00097     KalmanVertexFitter kvf(true);
00098     TransientVertex tv = kvf.vertex(transientTracks);
00099 
00100     if(tv.isValid()){
00101       GlobalVector tauDir(tau->px(),
00102           tau->py(),
00103           tau->pz());
00104       Vertex secVer = tv;
00105       // We have to un-const the PV for some reason
00106       reco::Vertex primaryVertexNonConst = *primaryVertex;
00107       flightPathSignificance = vertexSignificance(primaryVertexNonConst,secVer,tauDir);
00108     }
00109   }
00110   return flightPathSignificance;
00111 }
00112 
00113 double PFRecoTauDiscriminationByFlightPathSignificance::vertexSignificance(
00114     reco::Vertex& pv, Vertex& sv,GlobalVector& direction){
00115   return SecondaryVertex::computeDist3d(pv,sv,direction,withPVError).significance();
00116 }
00117 
00118 DEFINE_FWK_MODULE(PFRecoTauDiscriminationByFlightPathSignificance);
00119