00001 /* 00002 * ============================================================================= 00003 * Filename: RecoTauImpactParameterSignificancePlugin.cc 00004 * 00005 * Description: Add the IP significance of the lead track w.r.t to the PV. 00006 * to a PFTau. 00007 * Created: 10/31/2010 13:32:14 00008 * 00009 * Authors: Evan K. Friis (UC Davis), evan.klose.friis@cern.ch, 00010 * Simone Gennai, Ludovic Houchu 00011 * 00012 * ============================================================================= 00013 */ 00014 00015 #include "RecoTauTag/RecoTau/interface/RecoTauBuilderPlugins.h" 00016 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h" 00017 #include "DataFormats/VertexReco/interface/Vertex.h" 00018 #include "DataFormats/VertexReco/interface/VertexFwd.h" 00019 00020 #include "FWCore/Framework/interface/ESHandle.h" 00021 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h" 00022 00023 #include "TrackingTools/IPTools/interface/IPTools.h" 00024 #include "TrackingTools/Records/interface/TransientTrackRecord.h" 00025 #include "TrackingTools/TransientTrack/interface/TransientTrack.h" 00026 00027 namespace reco { namespace tau { 00028 00029 class RecoTauImpactParameterSignificancePlugin : public RecoTauModifierPlugin { 00030 public: 00031 explicit RecoTauImpactParameterSignificancePlugin( 00032 const edm::ParameterSet& pset); 00033 virtual ~RecoTauImpactParameterSignificancePlugin() {} 00034 void operator()(PFTau& tau) const; 00035 virtual void beginEvent(); 00036 private: 00037 edm::InputTag pvSrc_; 00038 const TransientTrackBuilder *builder_; 00039 const reco::Vertex* pv_; 00040 }; 00041 00042 RecoTauImpactParameterSignificancePlugin 00043 ::RecoTauImpactParameterSignificancePlugin(const edm::ParameterSet& pset) 00044 :RecoTauModifierPlugin(pset) { 00045 pvSrc_ = pset.getParameter<edm::InputTag>("pvSrc"); 00046 } 00047 00048 void RecoTauImpactParameterSignificancePlugin::beginEvent() { 00049 // Get primary vertex 00050 edm::Handle<reco::VertexCollection> pvs; 00051 evt()->getByLabel(pvSrc_, pvs); 00052 pv_ = &((*pvs)[0]); 00053 // Get tranisent track builder. 00054 edm::ESHandle<TransientTrackBuilder> myTransientTrackBuilder; 00055 evtSetup()->get<TransientTrackRecord>().get("TransientTrackBuilder", 00056 myTransientTrackBuilder); 00057 builder_= myTransientTrackBuilder.product(); 00058 } 00059 00060 void RecoTauImpactParameterSignificancePlugin::operator()(PFTau& tau) const { 00061 // Get the transient lead track 00062 if (tau.leadPFChargedHadrCand().isNonnull()) { 00063 TrackRef leadTrack = tau.leadPFChargedHadrCand()->trackRef(); 00064 if (leadTrack.isNonnull()) { 00065 const TransientTrack track = builder_->build(leadTrack); 00066 GlobalVector direction(tau.jetRef()->px(), tau.jetRef()->py(), 00067 tau.jetRef()->pz()); 00068 // Compute the significance 00069 std::pair<bool,Measurement1D> ipsig = 00070 //IPTools::signedTransverseImpactParameter(track, direction, *pv_); 00071 IPTools::signedImpactParameter3D(track, direction, *pv_); 00072 if (ipsig.first) 00073 tau.setleadPFChargedHadrCandsignedSipt(ipsig.second.significance()); 00074 } 00075 } 00076 } 00077 00078 }} // end namespace reco::tau 00079 #include "FWCore/Framework/interface/MakerMacros.h" 00080 DEFINE_EDM_PLUGIN(RecoTauModifierPluginFactory, 00081 reco::tau::RecoTauImpactParameterSignificancePlugin, 00082 "RecoTauImpactParameterSignificancePlugin");